xref: /petsc/src/ts/impls/explicit/rk/mrk.c (revision 474dd7730651a168e39c2c0e6a0673d505ddd4a4)
1*474dd773SHong Zhang /*
2*474dd773SHong Zhang   Code for time stepping with the multi-rate Runge-Kutta method
3*474dd773SHong Zhang 
4*474dd773SHong Zhang   Notes:
5*474dd773SHong Zhang   1) The general system is written as
6*474dd773SHong Zhang      Udot = F(t,U) for the nonsplit version of multi-rate RK,
7*474dd773SHong Zhang      user should give the indexes for both slow and fast components;
8*474dd773SHong Zhang   2) The general system is written as
9*474dd773SHong Zhang      Usdot = Fs(t,Us,Uf)
10*474dd773SHong Zhang      Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits,
11*474dd773SHong Zhang      user should partioned RHS by themselves and also provide the indexes for both slow and fast components.
12*474dd773SHong Zhang */
13*474dd773SHong Zhang 
14*474dd773SHong Zhang #include <petsc/private/tsimpl.h>
15*474dd773SHong Zhang #include <petscdm.h>
16*474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h>
17*474dd773SHong Zhang 
18*474dd773SHong Zhang static TSRKType TSMRKDefault = TSRK2A;
19*474dd773SHong Zhang 
20*474dd773SHong Zhang static PetscErrorCode TSSetUp_MRKSPLIT(TS ts)
21*474dd773SHong Zhang {
22*474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
23*474dd773SHong Zhang   RKTableau      tab = rk->tableau;
24*474dd773SHong Zhang   DM             dm,subdm,newdm;
25*474dd773SHong Zhang   PetscErrorCode ierr;
26*474dd773SHong Zhang 
27*474dd773SHong Zhang   PetscFunctionBegin;
28*474dd773SHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
29*474dd773SHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
30*474dd773SHong 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");
31*474dd773SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
32*474dd773SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
33*474dd773SHong 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");
34*474dd773SHong Zhang 
35*474dd773SHong Zhang   /* Only copy */
36*474dd773SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
37*474dd773SHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
38*474dd773SHong Zhang   ierr = TSGetDM(rk->subts_fast,&subdm);CHKERRQ(ierr);
39*474dd773SHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
40*474dd773SHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
41*474dd773SHong Zhang   ierr = TSSetDM(rk->subts_fast,newdm);CHKERRQ(ierr);
42*474dd773SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
43*474dd773SHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
44*474dd773SHong Zhang   ierr = TSGetDM(rk->subts_slow,&subdm);CHKERRQ(ierr);
45*474dd773SHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
46*474dd773SHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
47*474dd773SHong Zhang   ierr = TSSetDM(rk->subts_slow,newdm);CHKERRQ(ierr);
48*474dd773SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
49*474dd773SHong Zhang 
50*474dd773SHong Zhang   ierr = PetscMalloc1(tab->s,&rk->YdotRHS_fast);CHKERRQ(ierr);
51*474dd773SHong Zhang   ierr = PetscMalloc1(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
52*474dd773SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
53*474dd773SHong Zhang   PetscFunctionReturn(0);
54*474dd773SHong Zhang }
55*474dd773SHong Zhang 
56*474dd773SHong Zhang static PetscErrorCode TSReset_MRKSPLIT(TS ts)
57*474dd773SHong Zhang {
58*474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
59*474dd773SHong Zhang   PetscErrorCode ierr;
60*474dd773SHong Zhang 
61*474dd773SHong Zhang   PetscFunctionBegin;
62*474dd773SHong Zhang   ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr);
63*474dd773SHong Zhang   ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr);
64*474dd773SHong Zhang   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
65*474dd773SHong Zhang   PetscFunctionReturn(0);
66*474dd773SHong Zhang }
67*474dd773SHong Zhang 
68*474dd773SHong Zhang static PetscErrorCode TSSetUp_MRKNONSPLIT(TS ts)
69*474dd773SHong Zhang {
70*474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
71*474dd773SHong Zhang   RKTableau       tab  = rk->tableau;
72*474dd773SHong Zhang   PetscErrorCode ierr;
73*474dd773SHong Zhang 
74*474dd773SHong Zhang   PetscFunctionBegin;
75*474dd773SHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
76*474dd773SHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
77*474dd773SHong 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");
78*474dd773SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
79*474dd773SHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
80*474dd773SHong Zhang   PetscFunctionReturn(0);
81*474dd773SHong Zhang }
82*474dd773SHong Zhang 
83*474dd773SHong Zhang static PetscErrorCode TSReset_MRKNONSPLIT(TS ts)
84*474dd773SHong Zhang {
85*474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
86*474dd773SHong Zhang   RKTableau      tab  = rk->tableau;
87*474dd773SHong Zhang   PetscErrorCode ierr;
88*474dd773SHong Zhang 
89*474dd773SHong Zhang   PetscFunctionBegin;
90*474dd773SHong Zhang   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
91*474dd773SHong Zhang   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
92*474dd773SHong Zhang   PetscFunctionReturn(0);
93*474dd773SHong Zhang }
94*474dd773SHong Zhang 
95*474dd773SHong Zhang static PetscErrorCode TSInterpolate_MRKNONSPLIT(TS ts,PetscReal itime,Vec X)
96*474dd773SHong Zhang {
97*474dd773SHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
98*474dd773SHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
99*474dd773SHong Zhang   PetscReal        h;
100*474dd773SHong Zhang   PetscReal        tt,t;
101*474dd773SHong Zhang   PetscScalar      *b;
102*474dd773SHong Zhang   const PetscReal  *B = rk->tableau->binterp;
103*474dd773SHong Zhang   PetscErrorCode   ierr;
104*474dd773SHong Zhang 
105*474dd773SHong Zhang   PetscFunctionBegin;
106*474dd773SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
107*474dd773SHong Zhang 
108*474dd773SHong Zhang   switch (rk->status) {
109*474dd773SHong Zhang     case TS_STEP_INCOMPLETE:
110*474dd773SHong Zhang     case TS_STEP_PENDING:
111*474dd773SHong Zhang       h = ts->time_step;
112*474dd773SHong Zhang       t = (itime - ts->ptime)/h;
113*474dd773SHong Zhang       break;
114*474dd773SHong Zhang     case TS_STEP_COMPLETE:
115*474dd773SHong Zhang       h = ts->ptime - ts->ptime_prev;
116*474dd773SHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
117*474dd773SHong Zhang       break;
118*474dd773SHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
119*474dd773SHong Zhang   }
120*474dd773SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
121*474dd773SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
122*474dd773SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
123*474dd773SHong Zhang     for (i=0; i<s; i++) {
124*474dd773SHong Zhang       b[i]  += h * B[i*p+j] * tt;
125*474dd773SHong Zhang     }
126*474dd773SHong Zhang   }
127*474dd773SHong Zhang   ierr = VecCopy(rk->X0,X);CHKERRQ(ierr);
128*474dd773SHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
129*474dd773SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
130*474dd773SHong Zhang   PetscFunctionReturn(0);
131*474dd773SHong Zhang }
132*474dd773SHong Zhang 
133*474dd773SHong Zhang static PetscErrorCode TSInterpolate_MRKSPLIT(TS ts,PetscReal itime,Vec X)
134*474dd773SHong Zhang {
135*474dd773SHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
136*474dd773SHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
137*474dd773SHong Zhang   Vec              Yslow;    /* vector holds the slow components in Y[0] */
138*474dd773SHong Zhang   PetscReal        h;
139*474dd773SHong Zhang   PetscReal        tt,t;
140*474dd773SHong Zhang   PetscScalar      *b;
141*474dd773SHong Zhang   const PetscReal  *B = rk->tableau->binterp;
142*474dd773SHong Zhang   PetscErrorCode   ierr;
143*474dd773SHong Zhang 
144*474dd773SHong Zhang   PetscFunctionBegin;
145*474dd773SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
146*474dd773SHong Zhang 
147*474dd773SHong Zhang   switch (rk->status) {
148*474dd773SHong Zhang     case TS_STEP_INCOMPLETE:
149*474dd773SHong Zhang     case TS_STEP_PENDING:
150*474dd773SHong Zhang       h = ts->time_step;
151*474dd773SHong Zhang       t = (itime - ts->ptime)/h;
152*474dd773SHong Zhang       break;
153*474dd773SHong Zhang     case TS_STEP_COMPLETE:
154*474dd773SHong Zhang       h = ts->ptime - ts->ptime_prev;
155*474dd773SHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
156*474dd773SHong Zhang       break;
157*474dd773SHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
158*474dd773SHong Zhang   }
159*474dd773SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
160*474dd773SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
161*474dd773SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
162*474dd773SHong Zhang     for (i=0; i<s; i++) {
163*474dd773SHong Zhang       b[i]  += h * B[i*p+j] * tt;
164*474dd773SHong Zhang     }
165*474dd773SHong Zhang   }
166*474dd773SHong Zhang   ierr = VecGetSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr);
167*474dd773SHong Zhang   ierr = VecCopy(Yslow,X);CHKERRQ(ierr);
168*474dd773SHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
169*474dd773SHong Zhang   ierr = VecRestoreSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr);
170*474dd773SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
171*474dd773SHong Zhang   PetscFunctionReturn(0);
172*474dd773SHong Zhang }
173*474dd773SHong Zhang 
174*474dd773SHong Zhang static PetscErrorCode TSStep_MRKNONSPLIT(TS ts)
175*474dd773SHong Zhang {
176*474dd773SHong Zhang   TS_RK             *rk = (TS_RK*)ts->data;
177*474dd773SHong Zhang   RKTableau         tab  = rk->tableau;
178*474dd773SHong Zhang   Vec               *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow;
179*474dd773SHong Zhang   Vec               stage_slow,sol_slow;   /* vectors store the slow components */
180*474dd773SHong Zhang   Vec               subvec_slow;           /* sub vector to store the slow components */
181*474dd773SHong Zhang   const PetscInt    s = tab->s;
182*474dd773SHong Zhang   const PetscReal   *A = tab->A,*c = tab->c;
183*474dd773SHong Zhang   PetscScalar       *w = rk->work;
184*474dd773SHong Zhang   PetscInt          i,j,k,dtratio = rk->dtratio;
185*474dd773SHong Zhang   PetscReal         next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
186*474dd773SHong Zhang   PetscErrorCode    ierr;
187*474dd773SHong Zhang 
188*474dd773SHong Zhang   PetscFunctionBegin;
189*474dd773SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
190*474dd773SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr);
191*474dd773SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr);
192*474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
193*474dd773SHong Zhang   for (i=0; i<s; i++) {
194*474dd773SHong Zhang     rk->stage_time = t + h*c[i];
195*474dd773SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
196*474dd773SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
197*474dd773SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
198*474dd773SHong Zhang     ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr);
199*474dd773SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
200*474dd773SHong Zhang     /* compute the stage RHS */
201*474dd773SHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
202*474dd773SHong Zhang   }
203*474dd773SHong Zhang   /* update the slow components in the solution */
204*474dd773SHong Zhang   rk->YdotRHS = YdotRHS_slow;
205*474dd773SHong Zhang   rk->dtratio = 1;
206*474dd773SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr);
207*474dd773SHong Zhang   rk->dtratio = dtratio;
208*474dd773SHong Zhang   rk->YdotRHS = YdotRHS;
209*474dd773SHong Zhang   for (k=0; k<rk->dtratio; k++) {
210*474dd773SHong Zhang     for (i=0; i<s; i++) {
211*474dd773SHong Zhang       rk->stage_time = t + h/rk->dtratio*c[i];
212*474dd773SHong Zhang       ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
213*474dd773SHong 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 */
214*474dd773SHong Zhang       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
215*474dd773SHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
216*474dd773SHong Zhang       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
217*474dd773SHong Zhang       ierr = TSInterpolate_MRKNONSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],stage_slow);CHKERRQ(ierr);
218*474dd773SHong Zhang       /* update the slow components in the stage value */
219*474dd773SHong Zhang       ierr = VecGetSubVector(stage_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr);
220*474dd773SHong Zhang       ierr = VecISCopy(Y[i],rk->is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr);
221*474dd773SHong Zhang       ierr = VecRestoreSubVector(stage_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr);
222*474dd773SHong Zhang       ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr);
223*474dd773SHong Zhang       /* compute the stage RHS */
224*474dd773SHong Zhang       ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
225*474dd773SHong Zhang     }
226*474dd773SHong Zhang     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
227*474dd773SHong Zhang   }
228*474dd773SHong Zhang   /* update the slow components in the solution */
229*474dd773SHong Zhang   ierr = VecGetSubVector(sol_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr);
230*474dd773SHong Zhang   ierr = VecISCopy(ts->vec_sol,rk->is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr);
231*474dd773SHong Zhang   ierr = VecRestoreSubVector(sol_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr);
232*474dd773SHong Zhang 
233*474dd773SHong Zhang   ts->ptime += ts->time_step;
234*474dd773SHong Zhang   ts->time_step = next_time_step;
235*474dd773SHong Zhang   rk->status = TS_STEP_COMPLETE;
236*474dd773SHong Zhang   /* free memory of work vectors */
237*474dd773SHong Zhang   ierr = VecDestroy(&stage_slow);CHKERRQ(ierr);
238*474dd773SHong Zhang   ierr = VecDestroy(&sol_slow);CHKERRQ(ierr);
239*474dd773SHong Zhang   PetscFunctionReturn(0);
240*474dd773SHong Zhang }
241*474dd773SHong Zhang 
242*474dd773SHong Zhang /*
243*474dd773SHong Zhang  This is for partitioned RHS multirate RK method
244*474dd773SHong Zhang  The step completion formula is
245*474dd773SHong Zhang 
246*474dd773SHong Zhang  x1 = x0 + h b^T YdotRHS
247*474dd773SHong Zhang 
248*474dd773SHong Zhang */
249*474dd773SHong Zhang static PetscErrorCode TSEvaluateStep_MRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
250*474dd773SHong Zhang {
251*474dd773SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
252*474dd773SHong Zhang   RKTableau       tab  = rk->tableau;
253*474dd773SHong Zhang   Vec             Xslow,Xfast;                  /* subvectors of X which store slow components and fast components respectively */
254*474dd773SHong Zhang   PetscScalar     *w = rk->work;
255*474dd773SHong Zhang   PetscReal       h = ts->time_step;
256*474dd773SHong Zhang   PetscInt        s = tab->s,j;
257*474dd773SHong Zhang   PetscErrorCode  ierr;
258*474dd773SHong Zhang 
259*474dd773SHong Zhang   PetscFunctionBegin;
260*474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
261*474dd773SHong Zhang   if (rk->slow) {
262*474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h*tab->b[j];
263*474dd773SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
264*474dd773SHong Zhang     ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr);
265*474dd773SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);;
266*474dd773SHong Zhang   } else {
267*474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
268*474dd773SHong Zhang     ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
269*474dd773SHong Zhang     ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr);
270*474dd773SHong Zhang     ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
271*474dd773SHong Zhang   }
272*474dd773SHong Zhang   PetscFunctionReturn(0);
273*474dd773SHong Zhang }
274*474dd773SHong Zhang 
275*474dd773SHong Zhang static PetscErrorCode TSStep_MRKSPLIT(TS ts)
276*474dd773SHong Zhang {
277*474dd773SHong Zhang   TS_RK             *rk = (TS_RK*)ts->data;
278*474dd773SHong Zhang   RKTableau         tab = rk->tableau;
279*474dd773SHong Zhang   Vec               *Y = rk->Y,*YdotRHS = rk->YdotRHS;
280*474dd773SHong Zhang   Vec               *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow;
281*474dd773SHong Zhang   Vec               Yslow,Yfast;                         /* subvectors store the stges of slow components and fast components respectively                           */
282*474dd773SHong Zhang   const PetscInt    s = tab->s;
283*474dd773SHong Zhang   const PetscReal   *A = tab->A,*c = tab->c;
284*474dd773SHong Zhang   PetscScalar       *w = rk->work;
285*474dd773SHong Zhang   PetscInt          i,j,k;
286*474dd773SHong Zhang   PetscReal         next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
287*474dd773SHong Zhang   PetscErrorCode    ierr;
288*474dd773SHong Zhang 
289*474dd773SHong Zhang   PetscFunctionBegin;
290*474dd773SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
291*474dd773SHong Zhang   for (i=0; i<s; i++) {
292*474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
293*474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
294*474dd773SHong Zhang   }
295*474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
296*474dd773SHong Zhang   /* propogate both slow and fast components using large time steps */
297*474dd773SHong Zhang   for (i=0; i<s; i++) {
298*474dd773SHong Zhang     rk->stage_time = t + h*c[i];
299*474dd773SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
300*474dd773SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
301*474dd773SHong Zhang     /*ierr = VecCopy(ts->vec_sol,Yc[i]);CHKERRQ(ierr);*/
302*474dd773SHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
303*474dd773SHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
304*474dd773SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
305*474dd773SHong Zhang     ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr);
306*474dd773SHong Zhang     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
307*474dd773SHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
308*474dd773SHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
309*474dd773SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
310*474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
311*474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
312*474dd773SHong Zhang   }
313*474dd773SHong Zhang   rk->slow = PETSC_TRUE;
314*474dd773SHong Zhang   ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
315*474dd773SHong Zhang   for (k=0; k<rk->dtratio; k++) {
316*474dd773SHong Zhang     /* propogate fast component using small time steps */
317*474dd773SHong Zhang     for (i=0; i<s; i++) {
318*474dd773SHong Zhang       rk->stage_time = t + h/rk->dtratio*c[i];
319*474dd773SHong Zhang       ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
320*474dd773SHong Zhang       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
321*474dd773SHong Zhang       /* stage value for fast components */
322*474dd773SHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
323*474dd773SHong Zhang       ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
324*474dd773SHong Zhang       ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
325*474dd773SHong Zhang       ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
326*474dd773SHong Zhang       /* stage value for slow components */
327*474dd773SHong Zhang       ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
328*474dd773SHong Zhang       ierr = TSInterpolate_MRKSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Yslow);CHKERRQ(ierr);
329*474dd773SHong Zhang       ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
330*474dd773SHong Zhang       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
331*474dd773SHong Zhang       /* compute the stage RHS for fast components */
332*474dd773SHong Zhang       ierr = TSComputeRHSFunction(rk->subts_fast,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
333*474dd773SHong Zhang     }
334*474dd773SHong Zhang     /* update the value of fast components when using fast time step */
335*474dd773SHong Zhang     rk->slow = PETSC_FALSE;
336*474dd773SHong Zhang     ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
337*474dd773SHong Zhang   }
338*474dd773SHong Zhang   for (i=0; i<s; i++) {
339*474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
340*474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
341*474dd773SHong Zhang   }
342*474dd773SHong Zhang   ts->ptime += ts->time_step;
343*474dd773SHong Zhang   ts->time_step = next_time_step;
344*474dd773SHong Zhang   rk->status = TS_STEP_COMPLETE;
345*474dd773SHong Zhang   PetscFunctionReturn(0);
346*474dd773SHong Zhang }
347*474dd773SHong Zhang 
348*474dd773SHong Zhang /*@C
349*474dd773SHong Zhang   TSRKSetMultirateType - Set the type of RK Multirate scheme
350*474dd773SHong Zhang 
351*474dd773SHong Zhang   Logically collective
352*474dd773SHong Zhang 
353*474dd773SHong Zhang   Input Parameter:
354*474dd773SHong Zhang +  ts - timestepping context
355*474dd773SHong Zhang -  mrktype - type of MRK-scheme
356*474dd773SHong Zhang 
357*474dd773SHong Zhang   Options Database:
358*474dd773SHong Zhang .   -ts_rk_multiarte_type - <none,nonsplit,split>
359*474dd773SHong Zhang 
360*474dd773SHong Zhang   Level: intermediate
361*474dd773SHong Zhang @*/
362*474dd773SHong Zhang PetscErrorCode TSRKSetMultirateType(TS ts, TSMRKType mrktype)
363*474dd773SHong Zhang {
364*474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
365*474dd773SHong Zhang   PetscErrorCode ierr;
366*474dd773SHong Zhang 
367*474dd773SHong Zhang   PetscFunctionBegin;
368*474dd773SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
369*474dd773SHong Zhang   switch(mrktype){
370*474dd773SHong Zhang     case TSMRKNONE:
371*474dd773SHong Zhang       break;
372*474dd773SHong Zhang     case TSMRKNONSPLIT:
373*474dd773SHong Zhang       ts->ops->step           = TSStep_MRKNONSPLIT;
374*474dd773SHong Zhang       ts->ops->interpolate    = TSInterpolate_MRKNONSPLIT;
375*474dd773SHong Zhang       rk->dtratio             = 2;
376*474dd773SHong Zhang       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
377*474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKNONSPLIT_C",TSSetUp_MRKNONSPLIT);CHKERRQ(ierr);
378*474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKNONSPLIT_C",TSReset_MRKNONSPLIT);CHKERRQ(ierr);
379*474dd773SHong Zhang       break;
380*474dd773SHong Zhang     case TSMRKSPLIT:
381*474dd773SHong Zhang       ts->ops->step           = TSStep_MRKSPLIT;
382*474dd773SHong Zhang       ts->ops->evaluatestep   = TSEvaluateStep_MRKSPLIT;
383*474dd773SHong Zhang       ts->ops->interpolate    = TSInterpolate_MRKSPLIT;
384*474dd773SHong Zhang       rk->dtratio             = 2;
385*474dd773SHong Zhang       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
386*474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKSPLIT_C",TSSetUp_MRKSPLIT);CHKERRQ(ierr);
387*474dd773SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKSPLIT_C",TSReset_MRKSPLIT);CHKERRQ(ierr);
388*474dd773SHong Zhang       break;
389*474dd773SHong Zhang     default :
390*474dd773SHong Zhang       SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mrktype);
391*474dd773SHong Zhang   }
392*474dd773SHong Zhang   PetscFunctionReturn(0);
393*474dd773SHong Zhang }
394*474dd773SHong Zhang 
395