xref: /petsc/src/ts/impls/explicit/rk/mrk.c (revision feb237ba1401b79cfa9a99e81307a265c4a68462)
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>
1721052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h>
18474dd773SHong Zhang 
1902555c94SHong Zhang static PetscErrorCode TSReset_RK_MultirateNonsplit(TS ts)
20e5bffa22SHong Zhang {
21e5bffa22SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
22e5bffa22SHong Zhang   RKTableau      tab = rk->tableau;
23e5bffa22SHong Zhang   PetscErrorCode ierr;
24e5bffa22SHong Zhang 
25e5bffa22SHong Zhang   PetscFunctionBegin;
26e5bffa22SHong Zhang   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
27e5bffa22SHong Zhang   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
28e5bffa22SHong Zhang   PetscFunctionReturn(0);
29e5bffa22SHong Zhang }
30e5bffa22SHong Zhang 
3102555c94SHong Zhang static PetscErrorCode TSInterpolate_RK_MultirateNonsplit(TS ts,PetscReal itime,Vec X)
32e5bffa22SHong Zhang {
33e5bffa22SHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
34e5bffa22SHong Zhang   PetscInt         s = rk->tableau->s,p = rk->tableau->p,i,j;
3563a6f1b4SHong Zhang   PetscReal        h = ts->time_step;
36e5bffa22SHong Zhang   PetscReal        tt,t;
37e5bffa22SHong Zhang   PetscScalar      *b;
38e5bffa22SHong Zhang   const PetscReal  *B = rk->tableau->binterp;
39e5bffa22SHong Zhang   PetscErrorCode   ierr;
40e5bffa22SHong Zhang 
41e5bffa22SHong Zhang   PetscFunctionBegin;
42e5bffa22SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
4363a6f1b4SHong Zhang   t = (itime - rk->ptime)/h;
44e5bffa22SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
45e5bffa22SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
46e5bffa22SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
47e5bffa22SHong Zhang     for (i=0; i<s; i++) {
48e5bffa22SHong Zhang       b[i]  += h * B[i*p+j] * tt;
49e5bffa22SHong Zhang     }
50e5bffa22SHong Zhang   }
51e5bffa22SHong Zhang   ierr = VecCopy(rk->X0,X);CHKERRQ(ierr);
52e5bffa22SHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
53e5bffa22SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
54e5bffa22SHong Zhang   PetscFunctionReturn(0);
55e5bffa22SHong Zhang }
56e5bffa22SHong Zhang 
5702555c94SHong Zhang static PetscErrorCode TSStepRefine_RK_MultirateNonsplit(TS ts)
5863a6f1b4SHong Zhang {
5963a6f1b4SHong Zhang   TS              previousts,subts;
6063a6f1b4SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
6163a6f1b4SHong Zhang   RKTableau       tab = rk->tableau;
6263a6f1b4SHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
6363a6f1b4SHong Zhang   Vec             vec_fast,subvec_fast;
6463a6f1b4SHong Zhang   const PetscInt  s = tab->s;
6563a6f1b4SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
6663a6f1b4SHong Zhang   PetscScalar     *w = rk->work;
6763a6f1b4SHong Zhang   PetscInt        i,j,k;
6863a6f1b4SHong Zhang   PetscReal       t = ts->ptime,h = ts->time_step;
6963a6f1b4SHong Zhang   PetscErrorCode  ierr;
7063a6f1b4SHong Zhang 
7163a6f1b4SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&vec_fast);CHKERRQ(ierr);
7263a6f1b4SHong Zhang   previousts = rk->subts_current;
7363a6f1b4SHong Zhang   ierr = TSRHSSplitGetSubTS(rk->subts_current,"fast",&subts);CHKERRQ(ierr);
7463a6f1b4SHong Zhang   ierr = TSRHSSplitGetSubTS(subts,"fast",&subts);CHKERRQ(ierr);
7563a6f1b4SHong Zhang   for (k=0; k<rk->dtratio; k++) {
7663a6f1b4SHong Zhang     for (i=0; i<s; i++) {
7702555c94SHong Zhang       ierr = TSInterpolate_RK_MultirateNonsplit(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
7863a6f1b4SHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
7963a6f1b4SHong 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 */
8063a6f1b4SHong Zhang       ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr);
8163a6f1b4SHong Zhang       ierr = VecMAXPY(vec_fast,i,w,YdotRHS);CHKERRQ(ierr);
8263a6f1b4SHong Zhang       /* update the fast components in the stage value */
8363a6f1b4SHong Zhang       ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
8463a6f1b4SHong Zhang       ierr = VecISCopy(Y[i],rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr);
8563a6f1b4SHong Zhang       ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
8663a6f1b4SHong Zhang       /* compute the stage RHS */
8763a6f1b4SHong Zhang       ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
8863a6f1b4SHong Zhang     }
8963a6f1b4SHong Zhang     ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr);
9063a6f1b4SHong Zhang     ierr = TSEvaluateStep(ts,tab->order,vec_fast,NULL);CHKERRQ(ierr);
9163a6f1b4SHong Zhang     /* update the fast components in the solution */
9263a6f1b4SHong Zhang     ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
9363a6f1b4SHong Zhang     ierr = VecISCopy(ts->vec_sol,rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr);
9463a6f1b4SHong Zhang     ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
9563a6f1b4SHong Zhang 
9663a6f1b4SHong Zhang     if (subts) {
9763a6f1b4SHong Zhang       Vec *YdotRHS_copy;
9863a6f1b4SHong Zhang       ierr = VecDuplicateVecs(ts->vec_sol,s,&YdotRHS_copy);CHKERRQ(ierr);
9963a6f1b4SHong Zhang       rk->subts_current = rk->subts_fast;
10063a6f1b4SHong Zhang       ts->ptime = t+k*h/rk->dtratio;
10163a6f1b4SHong Zhang       ts->time_step = h/rk->dtratio;
10263a6f1b4SHong Zhang       ierr = TSRHSSplitGetIS(rk->subts_current,"fast",&rk->is_fast);CHKERRQ(ierr);
10363a6f1b4SHong Zhang       for (i=0; i<s; i++) {
10463a6f1b4SHong Zhang         ierr = VecCopy(rk->YdotRHS_slow[i],YdotRHS_copy[i]);CHKERRQ(ierr);
10563a6f1b4SHong Zhang         ierr = VecCopy(YdotRHS[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr);
10663a6f1b4SHong Zhang       }
10763a6f1b4SHong Zhang 
10802555c94SHong Zhang       ierr = TSStepRefine_RK_MultirateNonsplit(ts);CHKERRQ(ierr);
10963a6f1b4SHong Zhang 
11063a6f1b4SHong Zhang       rk->subts_current = previousts;
11163a6f1b4SHong Zhang       ts->ptime = t;
11263a6f1b4SHong Zhang       ts->time_step = h;
11363a6f1b4SHong Zhang       ierr = TSRHSSplitGetIS(previousts,"fast",&rk->is_fast);CHKERRQ(ierr);
11463a6f1b4SHong Zhang       for (i=0; i<s; i++) {
11563a6f1b4SHong Zhang         ierr = VecCopy(YdotRHS_copy[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr);
11663a6f1b4SHong Zhang       }
11763a6f1b4SHong Zhang       ierr = VecDestroyVecs(s,&YdotRHS_copy);CHKERRQ(ierr);
11863a6f1b4SHong Zhang     }
11963a6f1b4SHong Zhang   }
12063a6f1b4SHong Zhang   ierr = VecDestroy(&vec_fast);CHKERRQ(ierr);
12163a6f1b4SHong Zhang   PetscFunctionReturn(0);
12263a6f1b4SHong Zhang }
12363a6f1b4SHong Zhang 
12402555c94SHong Zhang static PetscErrorCode TSStep_RK_MultirateNonsplit(TS ts)
125e5bffa22SHong Zhang {
126e5bffa22SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
127e5bffa22SHong Zhang   RKTableau       tab = rk->tableau;
128e5bffa22SHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow;
129e5bffa22SHong Zhang   Vec             stage_slow,sol_slow; /* vectors store the slow components */
130e5bffa22SHong Zhang   Vec             subvec_slow; /* sub vector to store the slow components */
131e5bffa22SHong Zhang   IS              is_slow = rk->is_slow;
132e5bffa22SHong Zhang   const PetscInt  s = tab->s;
133e5bffa22SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
134e5bffa22SHong Zhang   PetscScalar     *w = rk->work;
13563a6f1b4SHong Zhang   PetscInt        i,j,dtratio = rk->dtratio;
136e5bffa22SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
137e5bffa22SHong Zhang   PetscErrorCode  ierr;
138e5bffa22SHong Zhang 
139e5bffa22SHong Zhang   PetscFunctionBegin;
140e5bffa22SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
141e5bffa22SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr);
142e5bffa22SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr);
143e5bffa22SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
144e5bffa22SHong Zhang   for (i=0; i<s; i++) {
145e5bffa22SHong Zhang     rk->stage_time = t + h*c[i];
146e5bffa22SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
147e5bffa22SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
148e5bffa22SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
149e5bffa22SHong Zhang     ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr);
150e5bffa22SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
151e5bffa22SHong Zhang     /* compute the stage RHS */
152e5bffa22SHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
153e5bffa22SHong Zhang   }
154e5bffa22SHong Zhang   /* update the slow components in the solution */
155e5bffa22SHong Zhang   rk->YdotRHS = YdotRHS_slow;
156e5bffa22SHong Zhang   rk->dtratio = 1;
157e5bffa22SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr);
158e5bffa22SHong Zhang   rk->dtratio = dtratio;
159e5bffa22SHong Zhang   rk->YdotRHS = YdotRHS;
160e5bffa22SHong Zhang   /* update the slow components in the solution */
161e5bffa22SHong Zhang   ierr = VecGetSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
162e5bffa22SHong Zhang   ierr = VecISCopy(ts->vec_sol,is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr);
163e5bffa22SHong Zhang   ierr = VecRestoreSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
164e5bffa22SHong Zhang 
16563a6f1b4SHong Zhang   rk->subts_current = ts;
16663a6f1b4SHong Zhang   rk->ptime = t;
16763a6f1b4SHong Zhang   rk->time_step = h;
16802555c94SHong Zhang   ierr = TSStepRefine_RK_MultirateNonsplit(ts);CHKERRQ(ierr);
16963a6f1b4SHong Zhang 
170e5bffa22SHong Zhang   ts->ptime = t + ts->time_step;
171e5bffa22SHong Zhang   ts->time_step = next_time_step;
172e5bffa22SHong Zhang   rk->status = TS_STEP_COMPLETE;
17363a6f1b4SHong Zhang 
174e5bffa22SHong Zhang   /* free memory of work vectors */
175e5bffa22SHong Zhang   ierr = VecDestroy(&stage_slow);CHKERRQ(ierr);
176e5bffa22SHong Zhang   ierr = VecDestroy(&sol_slow);CHKERRQ(ierr);
177e5bffa22SHong Zhang   PetscFunctionReturn(0);
178e5bffa22SHong Zhang }
179e5bffa22SHong Zhang 
18002555c94SHong Zhang static PetscErrorCode TSSetUp_RK_MultirateNonsplit(TS ts)
1810fe4d17eSHong Zhang {
1820fe4d17eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1830fe4d17eSHong Zhang   RKTableau      tab = rk->tableau;
1840fe4d17eSHong Zhang   PetscErrorCode ierr;
1850fe4d17eSHong Zhang 
1860fe4d17eSHong Zhang   PetscFunctionBegin;
1870fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
1880fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
1890fe4d17eSHong 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");
1900fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
1910fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
1920fe4d17eSHong 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");
1930fe4d17eSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
1940fe4d17eSHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
1950fe4d17eSHong Zhang   rk->subts_current = rk->subts_fast;
1960fe4d17eSHong Zhang 
19702555c94SHong Zhang   ts->ops->step        = TSStep_RK_MultirateNonsplit;
19802555c94SHong Zhang   ts->ops->interpolate = TSInterpolate_RK_MultirateNonsplit;
1990fe4d17eSHong Zhang   PetscFunctionReturn(0);
2000fe4d17eSHong Zhang }
2010fe4d17eSHong Zhang 
20263a6f1b4SHong Zhang /*
20363a6f1b4SHong Zhang   Copy DM from tssrc to tsdest, while keeping the original DMTS and DMSNES in tsdest.
20463a6f1b4SHong Zhang */
20563a6f1b4SHong Zhang static PetscErrorCode TSCopyDM(TS tssrc,TS tsdest)
20663a6f1b4SHong Zhang {
20763a6f1b4SHong Zhang   DM             newdm,dmsrc,dmdest;
20863a6f1b4SHong Zhang   PetscErrorCode ierr;
20963a6f1b4SHong Zhang 
21063a6f1b4SHong Zhang   PetscFunctionBegin;
21163a6f1b4SHong Zhang   ierr = TSGetDM(tssrc,&dmsrc);CHKERRQ(ierr);
21263a6f1b4SHong Zhang   ierr = DMClone(dmsrc,&newdm);CHKERRQ(ierr);
21363a6f1b4SHong Zhang   ierr = TSGetDM(tsdest,&dmdest);CHKERRQ(ierr);
21463a6f1b4SHong Zhang   ierr = DMCopyDMTS(dmdest,newdm);CHKERRQ(ierr);
21563a6f1b4SHong Zhang   ierr = DMCopyDMSNES(dmdest,newdm);CHKERRQ(ierr);
21663a6f1b4SHong Zhang   ierr = TSSetDM(tsdest,newdm);CHKERRQ(ierr);
21763a6f1b4SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
21863a6f1b4SHong Zhang   PetscFunctionReturn(0);
21963a6f1b4SHong Zhang }
22063a6f1b4SHong Zhang 
22102555c94SHong Zhang static PetscErrorCode TSReset_RK_MultirateSplit(TS ts)
222474dd773SHong Zhang {
223474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
224474dd773SHong Zhang   PetscErrorCode ierr;
225474dd773SHong Zhang 
226474dd773SHong Zhang   PetscFunctionBegin;
227bb42530cSHong Zhang   if (rk->subts_slow) {
228bb42530cSHong Zhang     ierr = PetscFree(rk->subts_slow->data);CHKERRQ(ierr);
229bb42530cSHong Zhang     rk->subts_slow = NULL;
230bb42530cSHong Zhang   }
231bb42530cSHong Zhang   if (rk->subts_fast) {
23263a6f1b4SHong Zhang     ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr);
23363a6f1b4SHong Zhang     ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr);
234bb42530cSHong Zhang     ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
235bf0cca7dSHong Zhang     ierr = TSReset_RK_MultirateSplit(rk->subts_fast);CHKERRQ(ierr);
236bb42530cSHong Zhang     ierr = PetscFree(rk->subts_fast->data);CHKERRQ(ierr);
237bb42530cSHong Zhang     rk->subts_fast = NULL;
238bb42530cSHong Zhang   }
239474dd773SHong Zhang   PetscFunctionReturn(0);
240474dd773SHong Zhang }
241474dd773SHong Zhang 
24202555c94SHong Zhang static PetscErrorCode TSInterpolate_RK_MultirateSplit(TS ts,PetscReal itime,Vec X)
243474dd773SHong Zhang {
244474dd773SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
24563a6f1b4SHong Zhang   Vec             Xslow;
246474dd773SHong Zhang   PetscInt        s = rk->tableau->s,p = rk->tableau->p,i,j;
247474dd773SHong Zhang   PetscReal       h;
248474dd773SHong Zhang   PetscReal       tt,t;
249474dd773SHong Zhang   PetscScalar     *b;
250474dd773SHong Zhang   const PetscReal *B = rk->tableau->binterp;
251474dd773SHong Zhang   PetscErrorCode  ierr;
252474dd773SHong Zhang 
253474dd773SHong Zhang   PetscFunctionBegin;
254474dd773SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
255474dd773SHong Zhang 
256474dd773SHong Zhang   switch (rk->status) {
257474dd773SHong Zhang     case TS_STEP_INCOMPLETE:
258474dd773SHong Zhang     case TS_STEP_PENDING:
259474dd773SHong Zhang       h = ts->time_step;
260474dd773SHong Zhang       t = (itime - ts->ptime)/h;
261474dd773SHong Zhang       break;
262474dd773SHong Zhang     case TS_STEP_COMPLETE:
263474dd773SHong Zhang       h = ts->ptime - ts->ptime_prev;
264474dd773SHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
265474dd773SHong Zhang       break;
266474dd773SHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
267474dd773SHong Zhang   }
268474dd773SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
269474dd773SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
270474dd773SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
271474dd773SHong Zhang     for (i=0; i<s; i++) {
272474dd773SHong Zhang       b[i]  += h * B[i*p+j] * tt;
273474dd773SHong Zhang     }
274474dd773SHong Zhang   }
27563a6f1b4SHong Zhang   for (i=0; i<s; i++) {
27663a6f1b4SHong Zhang     ierr = VecGetSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr);
27763a6f1b4SHong Zhang   }
27863a6f1b4SHong Zhang   ierr = VecGetSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr);
27963a6f1b4SHong Zhang   ierr = VecISCopy(rk->X0,rk->is_slow,SCATTER_REVERSE,Xslow);CHKERRQ(ierr);
28063a6f1b4SHong Zhang   ierr = VecMAXPY(Xslow,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
28163a6f1b4SHong Zhang   ierr = VecRestoreSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr);
28263a6f1b4SHong Zhang   for (i=0; i<s; i++) {
28363a6f1b4SHong Zhang     ierr = VecRestoreSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr);
28463a6f1b4SHong Zhang   }
285474dd773SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
286474dd773SHong Zhang   PetscFunctionReturn(0);
287474dd773SHong Zhang }
288474dd773SHong Zhang 
289474dd773SHong Zhang /*
290474dd773SHong Zhang  This is for partitioned RHS multirate RK method
291474dd773SHong Zhang  The step completion formula is
292474dd773SHong Zhang 
293474dd773SHong Zhang  x1 = x0 + h b^T YdotRHS
294474dd773SHong Zhang 
295474dd773SHong Zhang */
29602555c94SHong Zhang static PetscErrorCode TSEvaluateStep_RK_MultirateSplit(TS ts,PetscInt order,Vec X,PetscBool *done)
297474dd773SHong Zhang {
298474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
299474dd773SHong Zhang   RKTableau      tab = rk->tableau;
300474dd773SHong Zhang   Vec            Xslow,Xfast;                  /* subvectors of X which store slow components and fast components respectively */
301474dd773SHong Zhang   PetscScalar    *w = rk->work;
302474dd773SHong Zhang   PetscReal      h = ts->time_step;
303474dd773SHong Zhang   PetscInt       s = tab->s,j;
304474dd773SHong Zhang   PetscErrorCode ierr;
305474dd773SHong Zhang 
306474dd773SHong Zhang   PetscFunctionBegin;
307474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
308474dd773SHong Zhang   if (rk->slow) {
309474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h*tab->b[j];
310474dd773SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
31163a6f1b4SHong Zhang     ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr);
312*feb237baSPierre Jolivet     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
313474dd773SHong Zhang   } else {
314474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
315474dd773SHong Zhang     ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
31663a6f1b4SHong Zhang     ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr);
317474dd773SHong Zhang     ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
318474dd773SHong Zhang   }
319474dd773SHong Zhang   PetscFunctionReturn(0);
320474dd773SHong Zhang }
321474dd773SHong Zhang 
32202555c94SHong Zhang static PetscErrorCode TSStepRefine_RK_MultirateSplit(TS ts)
32363a6f1b4SHong Zhang {
32463a6f1b4SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
32563a6f1b4SHong Zhang   TS              subts_fast = rk->subts_fast,currentlevelts;
32663a6f1b4SHong Zhang   TS_RK           *subrk_fast = (TS_RK*)subts_fast->data;
32763a6f1b4SHong Zhang   RKTableau       tab = rk->tableau;
32863a6f1b4SHong Zhang   Vec             *Y = rk->Y;
32963a6f1b4SHong Zhang   Vec             *YdotRHS = rk->YdotRHS,*YdotRHS_fast = rk->YdotRHS_fast;
33063a6f1b4SHong Zhang   Vec             Yfast,Xfast;
33163a6f1b4SHong Zhang   const PetscInt  s = tab->s;
33263a6f1b4SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
33363a6f1b4SHong Zhang   PetscScalar     *w = rk->work;
33463a6f1b4SHong Zhang   PetscInt        i,j,k;
33563a6f1b4SHong Zhang   PetscReal       t = ts->ptime,h = ts->time_step;
33663a6f1b4SHong Zhang   PetscErrorCode  ierr;
33763a6f1b4SHong Zhang 
33863a6f1b4SHong Zhang   for (k=0; k<rk->dtratio; k++) {
33963a6f1b4SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
34063a6f1b4SHong Zhang     for (i=0; i<s; i++) {
34163a6f1b4SHong Zhang       ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
34263a6f1b4SHong Zhang     }
34363a6f1b4SHong Zhang     /* propogate fast component using small time steps */
34463a6f1b4SHong Zhang     for (i=0; i<s; i++) {
34563a6f1b4SHong Zhang       /* stage value for slow components */
34602555c94SHong Zhang       ierr = TSInterpolate_RK_MultirateSplit(rk->ts_root,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
34763a6f1b4SHong Zhang       currentlevelts = rk->ts_root;
34863a6f1b4SHong Zhang       while (currentlevelts != ts) { /* all the slow parts need to be interpolated separated */
34963a6f1b4SHong Zhang         currentlevelts = ((TS_RK*)currentlevelts->data)->subts_fast;
35002555c94SHong Zhang         ierr = TSInterpolate_RK_MultirateSplit(currentlevelts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
35163a6f1b4SHong Zhang       }
35263a6f1b4SHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
35363a6f1b4SHong Zhang       subrk_fast->stage_time = t + h/rk->dtratio*c[i];
35463a6f1b4SHong Zhang       ierr = TSPreStage(subts_fast,subrk_fast->stage_time);CHKERRQ(ierr);
35563a6f1b4SHong Zhang       /* stage value for fast components */
35663a6f1b4SHong Zhang       ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
35763a6f1b4SHong Zhang       ierr = VecCopy(Xfast,Yfast);CHKERRQ(ierr);
35863a6f1b4SHong Zhang       ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
35963a6f1b4SHong Zhang       ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
36063a6f1b4SHong Zhang       ierr = TSPostStage(subts_fast,subrk_fast->stage_time,i,Y);CHKERRQ(ierr);
36163a6f1b4SHong Zhang       /* compute the stage RHS for fast components */
36263a6f1b4SHong Zhang       ierr = TSComputeRHSFunction(subts_fast,t+k*h*rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
36363a6f1b4SHong Zhang     }
36463a6f1b4SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
36563a6f1b4SHong Zhang     /* update the value of fast components using fast time step */
36663a6f1b4SHong Zhang     rk->slow = PETSC_FALSE;
36702555c94SHong Zhang     ierr = TSEvaluateStep_RK_MultirateSplit(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
36863a6f1b4SHong Zhang     for (i=0; i<s; i++) {
36963a6f1b4SHong Zhang       ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
37063a6f1b4SHong Zhang     }
37163a6f1b4SHong Zhang 
37263a6f1b4SHong Zhang     if (subrk_fast->subts_fast) {
37363a6f1b4SHong Zhang       subts_fast->ptime = t+k*h/rk->dtratio;
37463a6f1b4SHong Zhang       subts_fast->time_step = h/rk->dtratio;
37502555c94SHong Zhang       ierr = TSStepRefine_RK_MultirateSplit(subts_fast);CHKERRQ(ierr);
37663a6f1b4SHong Zhang     }
37763a6f1b4SHong Zhang     /* update the fast components of the solution */
37863a6f1b4SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
37963a6f1b4SHong Zhang     ierr = VecISCopy(rk->X0,rk->is_fast,SCATTER_FORWARD,Xfast);CHKERRQ(ierr);
38063a6f1b4SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
38163a6f1b4SHong Zhang   }
38263a6f1b4SHong Zhang   PetscFunctionReturn(0);
38363a6f1b4SHong Zhang }
38463a6f1b4SHong Zhang 
38502555c94SHong Zhang static PetscErrorCode TSStep_RK_MultirateSplit(TS ts)
386474dd773SHong Zhang {
387474dd773SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
388474dd773SHong Zhang   RKTableau       tab = rk->tableau;
389474dd773SHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
39063a6f1b4SHong Zhang   Vec             *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow;
391474dd773SHong Zhang   Vec             Yslow,Yfast; /* subvectors store the stges of slow components and fast components respectively                           */
392474dd773SHong Zhang   const PetscInt  s = tab->s;
393474dd773SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
394474dd773SHong Zhang   PetscScalar     *w = rk->work;
39563a6f1b4SHong Zhang   PetscInt        i,j;
396474dd773SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
397474dd773SHong Zhang   PetscErrorCode  ierr;
398474dd773SHong Zhang 
399474dd773SHong Zhang   PetscFunctionBegin;
400474dd773SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
401474dd773SHong Zhang   for (i=0; i<s; i++) {
402474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
403474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
404474dd773SHong Zhang   }
405474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
406474dd773SHong Zhang   /* propogate both slow and fast components using large time steps */
407474dd773SHong Zhang   for (i=0; i<s; i++) {
408474dd773SHong Zhang     rk->stage_time = t + h*c[i];
409474dd773SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
410474dd773SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
411474dd773SHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
412bb42530cSHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
413474dd773SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
414474dd773SHong Zhang     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
415bb42530cSHong Zhang     ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr);
416474dd773SHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
417bb42530cSHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
418474dd773SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
419474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
420474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
421474dd773SHong Zhang   }
422474dd773SHong Zhang   rk->slow = PETSC_TRUE;
42363a6f1b4SHong Zhang   /* update the slow components of the solution using slow time step */
42402555c94SHong Zhang   ierr = TSEvaluateStep_RK_MultirateSplit(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
42563a6f1b4SHong Zhang   /* YdotRHS will be used for interpolation during refinement */
426474dd773SHong Zhang   for (i=0; i<s; i++) {
427474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
428474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
429474dd773SHong Zhang   }
43063a6f1b4SHong Zhang 
43102555c94SHong Zhang   ierr = TSStepRefine_RK_MultirateSplit(ts);CHKERRQ(ierr);
43263a6f1b4SHong Zhang 
433bb42530cSHong Zhang   ts->ptime = t + ts->time_step;
434474dd773SHong Zhang   ts->time_step = next_time_step;
435474dd773SHong Zhang   rk->status = TS_STEP_COMPLETE;
436474dd773SHong Zhang   PetscFunctionReturn(0);
437474dd773SHong Zhang }
438474dd773SHong Zhang 
43902555c94SHong Zhang static PetscErrorCode TSSetUp_RK_MultirateSplit(TS ts)
4400fe4d17eSHong Zhang {
4410fe4d17eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data,*nextlevelrk,*currentlevelrk;
4420fe4d17eSHong Zhang   TS             nextlevelts;
4430fe4d17eSHong Zhang   Vec            X0;
4440fe4d17eSHong Zhang   PetscErrorCode ierr;
4450fe4d17eSHong Zhang 
4460fe4d17eSHong Zhang   PetscFunctionBegin;
4470fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
4480fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
4490fe4d17eSHong 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");
4500fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
4510fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
4520fe4d17eSHong 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");
4530fe4d17eSHong Zhang 
4540fe4d17eSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&X0);CHKERRQ(ierr);
4550fe4d17eSHong Zhang   /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */
4560fe4d17eSHong Zhang   currentlevelrk = rk;
4570fe4d17eSHong Zhang   while (currentlevelrk->subts_fast) {
4580fe4d17eSHong Zhang     ierr = PetscMalloc1(rk->tableau->s,&currentlevelrk->YdotRHS_fast);CHKERRQ(ierr);
4590fe4d17eSHong Zhang     ierr = PetscMalloc1(rk->tableau->s,&currentlevelrk->YdotRHS_slow);CHKERRQ(ierr);
4600fe4d17eSHong Zhang     ierr = PetscObjectReference((PetscObject)X0);CHKERRQ(ierr);
4610fe4d17eSHong Zhang     currentlevelrk->X0 = X0;
4620fe4d17eSHong Zhang     currentlevelrk->ts_root = ts;
4630fe4d17eSHong Zhang 
4640fe4d17eSHong Zhang     /* set up the ts for the slow part */
4650fe4d17eSHong Zhang     nextlevelts = currentlevelrk->subts_slow;
4660fe4d17eSHong Zhang     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
4670fe4d17eSHong Zhang     nextlevelrk->tableau = rk->tableau;
4680fe4d17eSHong Zhang     nextlevelrk->work = rk->work;
4690fe4d17eSHong Zhang     nextlevelrk->Y = rk->Y;
4700fe4d17eSHong Zhang     nextlevelrk->YdotRHS = rk->YdotRHS;
4710fe4d17eSHong Zhang     nextlevelts->data = (void*)nextlevelrk;
4720fe4d17eSHong Zhang     ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr);
4730fe4d17eSHong Zhang     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
4740fe4d17eSHong Zhang 
4750fe4d17eSHong Zhang     /* set up the ts for the fast part */
4760fe4d17eSHong Zhang     nextlevelts = currentlevelrk->subts_fast;
4770fe4d17eSHong Zhang     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
4780fe4d17eSHong Zhang     nextlevelrk->tableau = rk->tableau;
4790fe4d17eSHong Zhang     nextlevelrk->work = rk->work;
4800fe4d17eSHong Zhang     nextlevelrk->Y = rk->Y;
4810fe4d17eSHong Zhang     nextlevelrk->YdotRHS = rk->YdotRHS;
4820fe4d17eSHong Zhang     nextlevelrk->dtratio = rk->dtratio;
4830fe4d17eSHong Zhang     ierr = TSRHSSplitGetIS(nextlevelts,"slow",&nextlevelrk->is_slow);CHKERRQ(ierr);
4840fe4d17eSHong Zhang     ierr = TSRHSSplitGetSubTS(nextlevelts,"slow",&nextlevelrk->subts_slow);CHKERRQ(ierr);
4850fe4d17eSHong Zhang     ierr = TSRHSSplitGetIS(nextlevelts,"fast",&nextlevelrk->is_fast);CHKERRQ(ierr);
4860fe4d17eSHong Zhang     ierr = TSRHSSplitGetSubTS(nextlevelts,"fast",&nextlevelrk->subts_fast);CHKERRQ(ierr);
4870fe4d17eSHong Zhang     nextlevelts->data = (void*)nextlevelrk;
4880fe4d17eSHong Zhang     ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr);
4890fe4d17eSHong Zhang     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
4900fe4d17eSHong Zhang 
4910fe4d17eSHong Zhang     currentlevelrk = nextlevelrk;
4920fe4d17eSHong Zhang   }
4930fe4d17eSHong Zhang   ierr = VecDestroy(&X0);CHKERRQ(ierr);
4940fe4d17eSHong Zhang 
49502555c94SHong Zhang   ts->ops->step         = TSStep_RK_MultirateSplit;
49602555c94SHong Zhang   ts->ops->evaluatestep = TSEvaluateStep_RK_MultirateSplit;
49702555c94SHong Zhang   ts->ops->interpolate  = TSInterpolate_RK_MultirateSplit;
4980fe4d17eSHong Zhang   PetscFunctionReturn(0);
4990fe4d17eSHong Zhang }
5000fe4d17eSHong Zhang 
50121052c0cSHong Zhang PetscErrorCode TSRKSetMultirate_RK(TS ts,PetscBool use_multirate)
502474dd773SHong Zhang {
503474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
504474dd773SHong Zhang   PetscErrorCode ierr;
505474dd773SHong Zhang 
506474dd773SHong Zhang   PetscFunctionBegin;
507474dd773SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5080fe4d17eSHong Zhang   rk->use_multirate = use_multirate;
5090fe4d17eSHong Zhang   if (use_multirate) {
510474dd773SHong Zhang     rk->dtratio = 2;
51102555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateSplit_C",TSSetUp_RK_MultirateSplit);CHKERRQ(ierr);
51202555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateSplit_C",TSReset_RK_MultirateSplit);CHKERRQ(ierr);
51302555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateNonsplit_C",TSSetUp_RK_MultirateNonsplit);CHKERRQ(ierr);
51402555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateNonsplit_C",TSReset_RK_MultirateNonsplit);CHKERRQ(ierr);
5150fe4d17eSHong Zhang   } else {
5160fe4d17eSHong Zhang     rk->dtratio = 0;
51702555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateSplit_C",NULL);CHKERRQ(ierr);
51802555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateSplit_C",NULL);CHKERRQ(ierr);
51902555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateNonsplit_C",NULL);CHKERRQ(ierr);
52002555c94SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateNonsplit_C",NULL);CHKERRQ(ierr);
521474dd773SHong Zhang   }
522474dd773SHong Zhang   PetscFunctionReturn(0);
523474dd773SHong Zhang }
524474dd773SHong Zhang 
52521052c0cSHong Zhang PetscErrorCode TSRKGetMultirate_RK(TS ts,PetscBool *use_multirate)
5260fe4d17eSHong Zhang {
5270fe4d17eSHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
5280fe4d17eSHong Zhang 
5290fe4d17eSHong Zhang   PetscFunctionBegin;
5300fe4d17eSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5310fe4d17eSHong Zhang   *use_multirate = rk->use_multirate;
5320fe4d17eSHong Zhang   PetscFunctionReturn(0);
5330fe4d17eSHong Zhang }
534