xref: /petsc/src/ts/impls/explicit/rk/mrk.c (revision e5bffa22434ee5e6fe35b8a1526b6ced76bd6539)
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 
18474dd773SHong Zhang static TSRKType TSMRKDefault = TSRK2A;
19474dd773SHong Zhang 
20*e5bffa22SHong Zhang static PetscErrorCode TSSetUp_MRKNONSPLIT(TS ts)
21*e5bffa22SHong Zhang {
22*e5bffa22SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
23*e5bffa22SHong Zhang   RKTableau       tab  = rk->tableau;
24*e5bffa22SHong Zhang   PetscErrorCode ierr;
25*e5bffa22SHong Zhang 
26*e5bffa22SHong Zhang   PetscFunctionBegin;
27*e5bffa22SHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
28*e5bffa22SHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
29*e5bffa22SHong 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");
30*e5bffa22SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
31*e5bffa22SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
32*e5bffa22SHong 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");
33*e5bffa22SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
34*e5bffa22SHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
35*e5bffa22SHong Zhang   PetscFunctionReturn(0);
36*e5bffa22SHong Zhang }
37*e5bffa22SHong Zhang 
38*e5bffa22SHong Zhang static PetscErrorCode TSReset_MRKNONSPLIT(TS ts)
39*e5bffa22SHong Zhang {
40*e5bffa22SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
41*e5bffa22SHong Zhang   RKTableau      tab  = rk->tableau;
42*e5bffa22SHong Zhang   PetscErrorCode ierr;
43*e5bffa22SHong Zhang 
44*e5bffa22SHong Zhang   PetscFunctionBegin;
45*e5bffa22SHong Zhang   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
46*e5bffa22SHong Zhang   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
47*e5bffa22SHong Zhang   PetscFunctionReturn(0);
48*e5bffa22SHong Zhang }
49*e5bffa22SHong Zhang 
50*e5bffa22SHong Zhang static PetscErrorCode TSInterpolate_MRKNONSPLIT(TS ts,PetscReal itime,Vec X)
51*e5bffa22SHong Zhang {
52*e5bffa22SHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
53*e5bffa22SHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
54*e5bffa22SHong Zhang   PetscReal        h;
55*e5bffa22SHong Zhang   PetscReal        tt,t;
56*e5bffa22SHong Zhang   PetscScalar      *b;
57*e5bffa22SHong Zhang   const PetscReal  *B = rk->tableau->binterp;
58*e5bffa22SHong Zhang   PetscErrorCode   ierr;
59*e5bffa22SHong Zhang 
60*e5bffa22SHong Zhang   PetscFunctionBegin;
61*e5bffa22SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
62*e5bffa22SHong Zhang 
63*e5bffa22SHong Zhang   switch (rk->status) {
64*e5bffa22SHong Zhang     case TS_STEP_INCOMPLETE:
65*e5bffa22SHong Zhang     case TS_STEP_PENDING:
66*e5bffa22SHong Zhang       h = ts->time_step;
67*e5bffa22SHong Zhang       t = (itime - ts->ptime)/h;
68*e5bffa22SHong Zhang       break;
69*e5bffa22SHong Zhang     case TS_STEP_COMPLETE:
70*e5bffa22SHong Zhang       h = ts->ptime - ts->ptime_prev;
71*e5bffa22SHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
72*e5bffa22SHong Zhang       break;
73*e5bffa22SHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
74*e5bffa22SHong Zhang   }
75*e5bffa22SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
76*e5bffa22SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
77*e5bffa22SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
78*e5bffa22SHong Zhang     for (i=0; i<s; i++) {
79*e5bffa22SHong Zhang       b[i]  += h * B[i*p+j] * tt;
80*e5bffa22SHong Zhang     }
81*e5bffa22SHong Zhang   }
82*e5bffa22SHong Zhang   ierr = VecCopy(rk->X0,X);CHKERRQ(ierr);
83*e5bffa22SHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
84*e5bffa22SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
85*e5bffa22SHong Zhang   PetscFunctionReturn(0);
86*e5bffa22SHong Zhang }
87*e5bffa22SHong Zhang 
88*e5bffa22SHong Zhang static PetscErrorCode TSStep_MRKNONSPLIT(TS ts)
89*e5bffa22SHong Zhang {
90*e5bffa22SHong Zhang   TS              nextlevelts;
91*e5bffa22SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
92*e5bffa22SHong Zhang   RKTableau       tab  = rk->tableau;
93*e5bffa22SHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow;
94*e5bffa22SHong Zhang   Vec             stage_slow,sol_slow;   /* vectors store the slow components */
95*e5bffa22SHong Zhang   Vec             subvec_slow;           /* sub vector to store the slow components */
96*e5bffa22SHong Zhang   IS              is_slow = rk->is_slow;
97*e5bffa22SHong Zhang   const PetscInt  s = tab->s;
98*e5bffa22SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
99*e5bffa22SHong Zhang   PetscScalar     *w = rk->work;
100*e5bffa22SHong Zhang   PetscInt        i,j,k,dtratio = rk->dtratio;
101*e5bffa22SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
102*e5bffa22SHong Zhang   PetscErrorCode  ierr;
103*e5bffa22SHong Zhang 
104*e5bffa22SHong Zhang   PetscFunctionBegin;
105*e5bffa22SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
106*e5bffa22SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr);
107*e5bffa22SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr);
108*e5bffa22SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
109*e5bffa22SHong Zhang   for (i=0; i<s; i++) {
110*e5bffa22SHong Zhang     rk->stage_time = t + h*c[i];
111*e5bffa22SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
112*e5bffa22SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
113*e5bffa22SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
114*e5bffa22SHong Zhang     ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr);
115*e5bffa22SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
116*e5bffa22SHong Zhang     /* compute the stage RHS */
117*e5bffa22SHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
118*e5bffa22SHong Zhang   }
119*e5bffa22SHong Zhang   /* update the slow components in the solution */
120*e5bffa22SHong Zhang   rk->YdotRHS = YdotRHS_slow;
121*e5bffa22SHong Zhang   rk->dtratio = 1;
122*e5bffa22SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr);
123*e5bffa22SHong Zhang   rk->dtratio = dtratio;
124*e5bffa22SHong Zhang   rk->YdotRHS = YdotRHS;
125*e5bffa22SHong Zhang 
126*e5bffa22SHong Zhang   ierr = TSRHSSplitGetSubTS(rk->subts_fast,"fast",&nextlevelts);CHKERRQ(ierr);
127*e5bffa22SHong Zhang   if (nextlevelts) {
128*e5bffa22SHong Zhang     for (k=0; k<rk->dtratio; k++) {
129*e5bffa22SHong Zhang       ierr = TSRHSSplitGetIS(nextlevelts,"slow",&rk->is_slow);CHKERRQ(ierr);
130*e5bffa22SHong Zhang       ts->time_step = h/rk->dtratio;
131*e5bffa22SHong Zhang       ierr = TSStep_MRKNONSPLIT(ts);CHKERRQ(ierr);
132*e5bffa22SHong Zhang       ts->time_step = h;
133*e5bffa22SHong Zhang     }
134*e5bffa22SHong Zhang   } else {
135*e5bffa22SHong Zhang     for (k=0; k<rk->dtratio; k++) {
136*e5bffa22SHong Zhang       for (i=0; i<s; i++) {
137*e5bffa22SHong Zhang         /* update the fast components in the stage value, the slow components will be overwritten, so it is ok to have garbage in the slow components */
138*e5bffa22SHong Zhang         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
139*e5bffa22SHong Zhang         for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
140*e5bffa22SHong Zhang         ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
141*e5bffa22SHong Zhang         ierr = TSInterpolate_MRKNONSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],stage_slow);CHKERRQ(ierr);
142*e5bffa22SHong Zhang         /* update the slow components in the stage value */
143*e5bffa22SHong Zhang         ierr = VecGetSubVector(stage_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr);
144*e5bffa22SHong Zhang         ierr = VecISCopy(Y[i],rk->is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr);
145*e5bffa22SHong Zhang         ierr = VecRestoreSubVector(stage_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr);
146*e5bffa22SHong Zhang         /* compute the stage RHS */
147*e5bffa22SHong Zhang         ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
148*e5bffa22SHong Zhang       }
149*e5bffa22SHong Zhang       ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
150*e5bffa22SHong Zhang     }
151*e5bffa22SHong Zhang   }
152*e5bffa22SHong Zhang   /* update the slow components in the solution */
153*e5bffa22SHong Zhang   ierr = VecGetSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
154*e5bffa22SHong Zhang   ierr = VecISCopy(ts->vec_sol,is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr);
155*e5bffa22SHong Zhang   ierr = VecRestoreSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
156*e5bffa22SHong Zhang 
157*e5bffa22SHong Zhang   ts->ptime = t + ts->time_step;
158*e5bffa22SHong Zhang   ts->time_step = next_time_step;
159*e5bffa22SHong Zhang   rk->status = TS_STEP_COMPLETE;
160*e5bffa22SHong Zhang   /* free memory of work vectors */
161*e5bffa22SHong Zhang   ierr = VecDestroy(&stage_slow);CHKERRQ(ierr);
162*e5bffa22SHong Zhang   ierr = VecDestroy(&sol_slow);CHKERRQ(ierr);
163*e5bffa22SHong Zhang   PetscFunctionReturn(0);
164*e5bffa22SHong Zhang }
165*e5bffa22SHong Zhang 
166474dd773SHong Zhang static PetscErrorCode TSSetUp_MRKSPLIT(TS ts)
167474dd773SHong Zhang {
168bb42530cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data,*nextlevelrk,*currentlevelrk;
169bb42530cSHong Zhang   TS             nextlevelts;
170474dd773SHong Zhang   DM             dm,subdm,newdm;
171474dd773SHong Zhang   PetscErrorCode ierr;
172474dd773SHong Zhang 
173474dd773SHong Zhang   PetscFunctionBegin;
174474dd773SHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
175474dd773SHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
176474dd773SHong 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");
177474dd773SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
178474dd773SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
179474dd773SHong 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");
180474dd773SHong Zhang 
181474dd773SHong Zhang   /* Only copy */
182474dd773SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
183474dd773SHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
184474dd773SHong Zhang   ierr = TSGetDM(rk->subts_fast,&subdm);CHKERRQ(ierr);
185474dd773SHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
186474dd773SHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
187474dd773SHong Zhang   ierr = TSSetDM(rk->subts_fast,newdm);CHKERRQ(ierr);
188474dd773SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
189474dd773SHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
190474dd773SHong Zhang   ierr = TSGetDM(rk->subts_slow,&subdm);CHKERRQ(ierr);
191474dd773SHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
192474dd773SHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
193474dd773SHong Zhang   ierr = TSSetDM(rk->subts_slow,newdm);CHKERRQ(ierr);
194474dd773SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
195474dd773SHong Zhang 
196474dd773SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
197bb42530cSHong Zhang   /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */
198bb42530cSHong Zhang   currentlevelrk = rk;
199bb42530cSHong Zhang   while (currentlevelrk->subts_fast) {
200bb42530cSHong Zhang     /* set up the ts for the slow part */
201bb42530cSHong Zhang     nextlevelts = currentlevelrk->subts_slow;
202bb42530cSHong Zhang     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
203bb42530cSHong Zhang     nextlevelrk->tableau = rk->tableau;
204bb42530cSHong Zhang     nextlevelrk->work = rk->work;
205bb42530cSHong Zhang     nextlevelrk->Y = rk->Y;
206bb42530cSHong Zhang     nextlevelrk->X0 = rk->X0;
207bb42530cSHong Zhang     ierr = PetscMalloc1(rk->tableau->s,&nextlevelrk->YdotRHS);CHKERRQ(ierr);
208bb42530cSHong Zhang     nextlevelts->data = (void*)nextlevelrk;
209bb42530cSHong Zhang     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
210bb42530cSHong Zhang 
211bb42530cSHong Zhang     /* set up the ts for the fast part */
212bb42530cSHong Zhang     nextlevelts = currentlevelrk->subts_fast;
213bb42530cSHong Zhang     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
214bb42530cSHong Zhang     nextlevelrk->tableau = rk->tableau;
215bb42530cSHong Zhang     nextlevelrk->work = rk->work;
216bb42530cSHong Zhang     nextlevelrk->Y = rk->Y;
217bb42530cSHong Zhang     nextlevelrk->X0 = rk->X0;
218bb42530cSHong Zhang     nextlevelrk->dtratio = rk->dtratio;
219bb42530cSHong Zhang     ierr = PetscMalloc1(rk->tableau->s,&nextlevelrk->YdotRHS);CHKERRQ(ierr);
220bb42530cSHong Zhang     ierr = TSRHSSplitGetIS(nextlevelts,"slow",&nextlevelrk->is_slow);CHKERRQ(ierr);
221bb42530cSHong Zhang     ierr = TSRHSSplitGetSubTS(nextlevelts,"slow",&nextlevelrk->subts_slow);CHKERRQ(ierr);
222bb42530cSHong Zhang     ierr = TSRHSSplitGetIS(nextlevelts,"fast",&nextlevelrk->is_fast);CHKERRQ(ierr);
223bb42530cSHong Zhang     ierr = TSRHSSplitGetSubTS(nextlevelts,"fast",&nextlevelrk->subts_fast);CHKERRQ(ierr);
224bb42530cSHong Zhang     nextlevelts->data = (void*)nextlevelrk;
225bb42530cSHong Zhang     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
226bb42530cSHong Zhang 
227bb42530cSHong Zhang     currentlevelrk = nextlevelrk;
228bb42530cSHong Zhang   }
229474dd773SHong Zhang   PetscFunctionReturn(0);
230474dd773SHong Zhang }
231474dd773SHong Zhang 
232474dd773SHong Zhang static PetscErrorCode TSReset_MRKSPLIT(TS ts)
233474dd773SHong Zhang {
234474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
235474dd773SHong Zhang   PetscErrorCode ierr;
236474dd773SHong Zhang 
237474dd773SHong Zhang   PetscFunctionBegin;
238bb42530cSHong Zhang   if (rk->subts_slow) {
239bb42530cSHong Zhang     ierr = TSReset_MRKSPLIT(rk->subts_slow);CHKERRQ(ierr);
240bb42530cSHong Zhang     ierr = PetscFree(rk->subts_slow->data);CHKERRQ(ierr);
241bb42530cSHong Zhang     rk->subts_slow = NULL;
242bb42530cSHong Zhang   }
243bb42530cSHong Zhang   if (rk->subts_fast) {
244bb42530cSHong Zhang     ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
245bb42530cSHong Zhang     ierr = TSReset_MRKSPLIT(rk->subts_fast);
246bb42530cSHong Zhang     ierr = PetscFree(rk->subts_fast->data);CHKERRQ(ierr);
247bb42530cSHong Zhang     rk->subts_fast = NULL;
248bb42530cSHong Zhang   }
249474dd773SHong Zhang   ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr);
250474dd773SHong Zhang   ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr);
251bb42530cSHong Zhang   ierr = PetscFree(rk->YdotRHS);CHKERRQ(ierr);
252474dd773SHong Zhang   PetscFunctionReturn(0);
253474dd773SHong Zhang }
254474dd773SHong Zhang 
255474dd773SHong Zhang static PetscErrorCode TSInterpolate_MRKSPLIT(TS ts,PetscReal itime,Vec X)
256474dd773SHong Zhang {
257474dd773SHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
258474dd773SHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
259474dd773SHong Zhang   Vec              Yslow;    /* vector holds the slow components in Y[0] */
260474dd773SHong Zhang   PetscReal        h;
261474dd773SHong Zhang   PetscReal        tt,t;
262474dd773SHong Zhang   PetscScalar      *b;
263474dd773SHong Zhang   const PetscReal  *B = rk->tableau->binterp;
264474dd773SHong Zhang   PetscErrorCode   ierr;
265474dd773SHong Zhang 
266474dd773SHong Zhang   PetscFunctionBegin;
267474dd773SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
268474dd773SHong Zhang 
269474dd773SHong Zhang   switch (rk->status) {
270474dd773SHong Zhang     case TS_STEP_INCOMPLETE:
271474dd773SHong Zhang     case TS_STEP_PENDING:
272474dd773SHong Zhang       h = ts->time_step;
273474dd773SHong Zhang       t = (itime - ts->ptime)/h;
274474dd773SHong Zhang       break;
275474dd773SHong Zhang     case TS_STEP_COMPLETE:
276474dd773SHong Zhang       h = ts->ptime - ts->ptime_prev;
277474dd773SHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
278474dd773SHong Zhang       break;
279474dd773SHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
280474dd773SHong Zhang   }
281474dd773SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
282474dd773SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
283474dd773SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
284474dd773SHong Zhang     for (i=0; i<s; i++) {
285474dd773SHong Zhang       b[i]  += h * B[i*p+j] * tt;
286474dd773SHong Zhang     }
287474dd773SHong Zhang   }
288474dd773SHong Zhang   ierr = VecGetSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr);
289474dd773SHong Zhang   ierr = VecCopy(Yslow,X);CHKERRQ(ierr);
290bb42530cSHong Zhang   ierr = VecMAXPY(X,s,b,((TS_RK*)rk->subts_slow->data)->YdotRHS);CHKERRQ(ierr);
291474dd773SHong Zhang   ierr = VecRestoreSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr);
292474dd773SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
293474dd773SHong Zhang   PetscFunctionReturn(0);
294474dd773SHong Zhang }
295474dd773SHong Zhang 
296474dd773SHong Zhang /*
297474dd773SHong Zhang  This is for partitioned RHS multirate RK method
298474dd773SHong Zhang  The step completion formula is
299474dd773SHong Zhang 
300474dd773SHong Zhang  x1 = x0 + h b^T YdotRHS
301474dd773SHong Zhang 
302474dd773SHong Zhang */
303474dd773SHong Zhang static PetscErrorCode TSEvaluateStep_MRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
304474dd773SHong Zhang {
305474dd773SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
306474dd773SHong Zhang   RKTableau       tab  = rk->tableau;
307474dd773SHong Zhang   Vec             Xslow,Xfast;                  /* subvectors of X which store slow components and fast components respectively */
308474dd773SHong Zhang   PetscScalar     *w = rk->work;
309474dd773SHong Zhang   PetscReal       h = ts->time_step;
310474dd773SHong Zhang   PetscInt        s = tab->s,j;
311474dd773SHong Zhang   PetscErrorCode  ierr;
312474dd773SHong Zhang 
313474dd773SHong Zhang   PetscFunctionBegin;
314474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
315474dd773SHong Zhang   if (rk->slow) {
316474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h*tab->b[j];
317474dd773SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
318bb42530cSHong Zhang     ierr = VecMAXPY(Xslow,s,w,((TS_RK*)rk->subts_slow->data)->YdotRHS);CHKERRQ(ierr);
319474dd773SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);;
320474dd773SHong Zhang   } else {
321474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
322474dd773SHong Zhang     ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
323bb42530cSHong Zhang     ierr = VecMAXPY(Xfast,s,w,((TS_RK*)rk->subts_fast->data)->YdotRHS);CHKERRQ(ierr);
324474dd773SHong Zhang     ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
325474dd773SHong Zhang   }
326474dd773SHong Zhang   PetscFunctionReturn(0);
327474dd773SHong Zhang }
328474dd773SHong Zhang 
329474dd773SHong Zhang static PetscErrorCode TSStep_MRKSPLIT(TS ts)
330474dd773SHong Zhang {
331474dd773SHong Zhang   TS_RK             *rk = (TS_RK*)ts->data;
332474dd773SHong Zhang   RKTableau         tab = rk->tableau;
333474dd773SHong Zhang   Vec               *Y = rk->Y,*YdotRHS = rk->YdotRHS;
334bb42530cSHong Zhang   Vec               *YdotRHS_fast,*YdotRHS_slow;
335474dd773SHong Zhang   Vec               Yslow,Yfast;                         /* subvectors store the stges of slow components and fast components respectively                           */
336474dd773SHong Zhang   const PetscInt    s = tab->s;
337474dd773SHong Zhang   const PetscReal   *A = tab->A,*c = tab->c;
338474dd773SHong Zhang   PetscScalar       *w = rk->work;
339474dd773SHong Zhang   PetscInt          i,j,k;
340474dd773SHong Zhang   PetscReal         next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
341474dd773SHong Zhang   PetscErrorCode    ierr;
342474dd773SHong Zhang 
343474dd773SHong Zhang   PetscFunctionBegin;
344474dd773SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
345bb42530cSHong Zhang   YdotRHS_fast = ((TS_RK*)rk->subts_fast->data)->YdotRHS;
346bb42530cSHong Zhang   YdotRHS_slow = ((TS_RK*)rk->subts_slow->data)->YdotRHS;
347474dd773SHong Zhang   for (i=0; i<s; i++) {
348474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
349474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
350474dd773SHong Zhang   }
351474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
352474dd773SHong Zhang   /* propogate both slow and fast components using large time steps */
353474dd773SHong Zhang   for (i=0; i<s; i++) {
354474dd773SHong Zhang     rk->stage_time = t + h*c[i];
355474dd773SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
356474dd773SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
357474dd773SHong Zhang     /*ierr = VecCopy(ts->vec_sol,Yc[i]);CHKERRQ(ierr);*/
358474dd773SHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
359bb42530cSHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
360474dd773SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
361474dd773SHong Zhang     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
362bb42530cSHong Zhang     ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr);
363474dd773SHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
364bb42530cSHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
365474dd773SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
366474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
367474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
368474dd773SHong Zhang   }
369474dd773SHong Zhang   rk->slow = PETSC_TRUE;
370bb42530cSHong Zhang   /* update the value of slow components when using slow time step */
371474dd773SHong Zhang   ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
372bb42530cSHong Zhang 
373bb42530cSHong Zhang   if (((TS_RK*)rk->subts_fast->data)->subts_fast) {
374bb42530cSHong Zhang     rk->subts_fast->ptime = t;
375bb42530cSHong Zhang     rk->subts_fast->time_step = h/rk->dtratio;
376bb42530cSHong Zhang     for (k=0; k<rk->dtratio; k++) {
377bb42530cSHong Zhang       ierr = TSStep_MRKSPLIT(rk->subts_fast);CHKERRQ(ierr);
378bb42530cSHong Zhang     }
379bb42530cSHong Zhang   } else {
380474dd773SHong Zhang     for (k=0; k<rk->dtratio; k++) {
381474dd773SHong Zhang       /* propogate fast component using small time steps */
382474dd773SHong Zhang       for (i=0; i<s; i++) {
383474dd773SHong Zhang         rk->stage_time = t + h/rk->dtratio*c[i];
384474dd773SHong Zhang         ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
385474dd773SHong Zhang         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
386474dd773SHong Zhang         /* stage value for fast components */
387474dd773SHong Zhang         for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
388474dd773SHong Zhang         ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
389474dd773SHong Zhang         ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
390474dd773SHong Zhang         ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
391474dd773SHong Zhang         /* stage value for slow components */
392474dd773SHong Zhang         ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
393474dd773SHong Zhang         ierr = TSInterpolate_MRKSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Yslow);CHKERRQ(ierr);
394474dd773SHong Zhang         ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
395474dd773SHong Zhang         ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr);
396474dd773SHong Zhang         /* compute the stage RHS for fast components */
397bb42530cSHong Zhang         ierr = TSComputeRHSFunction(rk->subts_fast,t+k*h*rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
398474dd773SHong Zhang       }
399474dd773SHong Zhang       /* update the value of fast components when using fast time step */
400474dd773SHong Zhang       rk->slow = PETSC_FALSE;
401474dd773SHong Zhang       ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
402474dd773SHong Zhang     }
403bb42530cSHong Zhang   }
404474dd773SHong Zhang   for (i=0; i<s; i++) {
405474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
406474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
407474dd773SHong Zhang   }
408bb42530cSHong Zhang   ts->ptime = t + ts->time_step;
409474dd773SHong Zhang   ts->time_step = next_time_step;
410474dd773SHong Zhang   rk->status = TS_STEP_COMPLETE;
411474dd773SHong Zhang   PetscFunctionReturn(0);
412474dd773SHong Zhang }
413474dd773SHong Zhang 
414474dd773SHong Zhang /*@C
415474dd773SHong Zhang   TSRKSetMultirateType - Set the type of RK Multirate scheme
416474dd773SHong Zhang 
417474dd773SHong Zhang   Logically collective
418474dd773SHong Zhang 
419474dd773SHong Zhang   Input Parameter:
420474dd773SHong Zhang +  ts - timestepping context
421474dd773SHong Zhang -  mrktype - type of MRK-scheme
422474dd773SHong Zhang 
423474dd773SHong Zhang   Options Database:
424474dd773SHong Zhang .   -ts_rk_multiarte_type - <none,nonsplit,split>
425474dd773SHong Zhang 
426474dd773SHong Zhang   Level: intermediate
427474dd773SHong Zhang @*/
428474dd773SHong Zhang PetscErrorCode TSRKSetMultirateType(TS ts, TSMRKType mrktype)
429474dd773SHong Zhang {
430474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
431474dd773SHong Zhang   PetscErrorCode ierr;
432474dd773SHong Zhang 
433474dd773SHong Zhang   PetscFunctionBegin;
434474dd773SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
435474dd773SHong Zhang   switch(mrktype){
436474dd773SHong Zhang     case TSMRKNONE:
437474dd773SHong Zhang       break;
438474dd773SHong Zhang     case TSMRKNONSPLIT:
439474dd773SHong Zhang       ts->ops->step           = TSStep_MRKNONSPLIT;
440474dd773SHong Zhang       ts->ops->interpolate    = TSInterpolate_MRKNONSPLIT;
441474dd773SHong Zhang       rk->dtratio             = 2;
442474dd773SHong Zhang       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
443474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKNONSPLIT_C",TSSetUp_MRKNONSPLIT);CHKERRQ(ierr);
444474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKNONSPLIT_C",TSReset_MRKNONSPLIT);CHKERRQ(ierr);
445474dd773SHong Zhang       break;
446474dd773SHong Zhang     case TSMRKSPLIT:
447474dd773SHong Zhang       ts->ops->step           = TSStep_MRKSPLIT;
448474dd773SHong Zhang       ts->ops->evaluatestep   = TSEvaluateStep_MRKSPLIT;
449474dd773SHong Zhang       ts->ops->interpolate    = TSInterpolate_MRKSPLIT;
450474dd773SHong Zhang       rk->dtratio             = 2;
451474dd773SHong Zhang       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
452474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKSPLIT_C",TSSetUp_MRKSPLIT);CHKERRQ(ierr);
453474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKSPLIT_C",TSReset_MRKSPLIT);CHKERRQ(ierr);
454474dd773SHong Zhang       break;
455474dd773SHong Zhang     default :
456474dd773SHong Zhang       SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mrktype);
457474dd773SHong Zhang   }
458474dd773SHong Zhang   PetscFunctionReturn(0);
459474dd773SHong Zhang }
460474dd773SHong Zhang 
461