xref: /petsc/src/ts/impls/arkimex/fsarkimex.c (revision 3a2a065f407e7da3bb6fc7eaf619dbab56034fde)
1*3a2a065fSHong Zhang #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
2*3a2a065fSHong Zhang #include <petscdm.h>
3*3a2a065fSHong Zhang #include <../src/ts/impls/arkimex/arkimex.h>
4*3a2a065fSHong Zhang #include <../src/ts/impls/arkimex/fsarkimex.h>
5*3a2a065fSHong Zhang 
6*3a2a065fSHong Zhang static PetscErrorCode TSARKIMEXSetSplits(TS ts)
7*3a2a065fSHong Zhang {
8*3a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
9*3a2a065fSHong Zhang   DM          dm, subdm, newdm;
10*3a2a065fSHong Zhang 
11*3a2a065fSHong Zhang   PetscFunctionBegin;
12*3a2a065fSHong Zhang   PetscCall(TSRHSSplitGetSubTS(ts, "slow", &ark->subts_slow));
13*3a2a065fSHong Zhang   PetscCall(TSRHSSplitGetSubTS(ts, "fast", &ark->subts_fast));
14*3a2a065fSHong Zhang   /* Only copy the DM */
15*3a2a065fSHong Zhang   PetscCall(TSGetDM(ts, &dm));
16*3a2a065fSHong Zhang   if (ark->subts_slow) {
17*3a2a065fSHong Zhang     PetscCall(DMClone(dm, &newdm));
18*3a2a065fSHong Zhang     PetscCall(TSGetDM(ark->subts_slow, &subdm));
19*3a2a065fSHong Zhang     PetscCall(DMCopyDMTS(subdm, newdm));
20*3a2a065fSHong Zhang     PetscCall(TSSetDM(ark->subts_slow, newdm));
21*3a2a065fSHong Zhang     PetscCall(DMDestroy(&newdm));
22*3a2a065fSHong Zhang   }
23*3a2a065fSHong Zhang   if (ark->subts_fast) {
24*3a2a065fSHong Zhang     PetscCall(DMClone(dm, &newdm));
25*3a2a065fSHong Zhang     PetscCall(TSGetDM(ark->subts_fast, &subdm));
26*3a2a065fSHong Zhang     PetscCall(DMCopyDMTS(subdm, newdm));
27*3a2a065fSHong Zhang     PetscCall(TSSetDM(ark->subts_fast, newdm));
28*3a2a065fSHong Zhang     PetscCall(DMDestroy(&newdm));
29*3a2a065fSHong Zhang   }
30*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
31*3a2a065fSHong Zhang }
32*3a2a065fSHong Zhang 
33*3a2a065fSHong Zhang static PetscErrorCode TSExtrapolate_ARKIMEX_FastSlowSplit(TS ts, PetscReal c, Vec X)
34*3a2a065fSHong Zhang {
35*3a2a065fSHong Zhang   TS_ARKIMEX      *ark = (TS_ARKIMEX *)ts->data;
36*3a2a065fSHong Zhang   ARKTableau       tab = ark->tableau;
37*3a2a065fSHong Zhang   PetscInt         s = tab->s, pinterp = tab->pinterp, i, j;
38*3a2a065fSHong Zhang   PetscReal        h, h_prev, t, tt;
39*3a2a065fSHong Zhang   PetscScalar     *bt = ark->work, *b = ark->work + s;
40*3a2a065fSHong Zhang   const PetscReal *Bt = tab->binterpt, *B = tab->binterp;
41*3a2a065fSHong Zhang   PetscBool        fasthasE;
42*3a2a065fSHong Zhang 
43*3a2a065fSHong Zhang   PetscFunctionBegin;
44*3a2a065fSHong Zhang   PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name);
45*3a2a065fSHong Zhang   h      = ts->time_step;
46*3a2a065fSHong Zhang   h_prev = ts->ptime - ts->ptime_prev;
47*3a2a065fSHong Zhang   t      = 1 + h / h_prev * c;
48*3a2a065fSHong Zhang   for (i = 0; i < s; i++) bt[i] = b[i] = 0;
49*3a2a065fSHong Zhang   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
50*3a2a065fSHong Zhang     for (i = 0; i < s; i++) {
51*3a2a065fSHong Zhang       bt[i] += h * Bt[i * pinterp + j] * tt;
52*3a2a065fSHong Zhang       b[i] += h * B[i * pinterp + j] * tt;
53*3a2a065fSHong Zhang     }
54*3a2a065fSHong Zhang   }
55*3a2a065fSHong Zhang   PetscCheck(ark->Y_prev, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Stages from previous step have not been stored");
56*3a2a065fSHong Zhang   PetscCall(VecCopy(ark->Y_prev[0], X));
57*3a2a065fSHong Zhang   PetscCall(VecMAXPY(X, s, bt, ark->YdotI_prev));
58*3a2a065fSHong Zhang   PetscCall(TSHasRHSFunction(ark->subts_fast, &fasthasE));
59*3a2a065fSHong Zhang   if (fasthasE) PetscCall(VecMAXPY(X, s, b, ark->YdotRHS_prev));
60*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
61*3a2a065fSHong Zhang }
62*3a2a065fSHong Zhang 
63*3a2a065fSHong Zhang /*
64*3a2a065fSHong Zhang  The step completion formula is
65*3a2a065fSHong Zhang 
66*3a2a065fSHong Zhang  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
67*3a2a065fSHong Zhang 
68*3a2a065fSHong Zhang  This function can be called before or after ts->vec_sol has been updated.
69*3a2a065fSHong Zhang  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
70*3a2a065fSHong Zhang  We can write
71*3a2a065fSHong Zhang 
72*3a2a065fSHong Zhang  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
73*3a2a065fSHong Zhang      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
74*3a2a065fSHong Zhang      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
75*3a2a065fSHong Zhang 
76*3a2a065fSHong Zhang  so we can evaluate the method with different order even after the step has been optimistically completed.
77*3a2a065fSHong Zhang */
78*3a2a065fSHong Zhang static PetscErrorCode TSEvaluateStep_ARKIMEX_FastSlowSplit(TS ts, PetscInt order, Vec X, PetscBool *done)
79*3a2a065fSHong Zhang {
80*3a2a065fSHong Zhang   TS_ARKIMEX  *ark = (TS_ARKIMEX *)ts->data;
81*3a2a065fSHong Zhang   ARKTableau   tab = ark->tableau;
82*3a2a065fSHong Zhang   Vec          Xfast, Xslow;
83*3a2a065fSHong Zhang   PetscScalar *w = ark->work;
84*3a2a065fSHong Zhang   PetscReal    h;
85*3a2a065fSHong Zhang   PetscInt     s = tab->s, j;
86*3a2a065fSHong Zhang   PetscBool    fasthasE;
87*3a2a065fSHong Zhang 
88*3a2a065fSHong Zhang   PetscFunctionBegin;
89*3a2a065fSHong Zhang   switch (ark->status) {
90*3a2a065fSHong Zhang   case TS_STEP_INCOMPLETE:
91*3a2a065fSHong Zhang   case TS_STEP_PENDING:
92*3a2a065fSHong Zhang     h = ts->time_step;
93*3a2a065fSHong Zhang     break;
94*3a2a065fSHong Zhang   case TS_STEP_COMPLETE:
95*3a2a065fSHong Zhang     h = ts->ptime - ts->ptime_prev;
96*3a2a065fSHong Zhang     break;
97*3a2a065fSHong Zhang   default:
98*3a2a065fSHong Zhang     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
99*3a2a065fSHong Zhang   }
100*3a2a065fSHong Zhang   if (ark->is_fast) PetscCall(TSHasRHSFunction(ark->subts_fast, &fasthasE));
101*3a2a065fSHong Zhang   if (order == tab->order) {
102*3a2a065fSHong Zhang     if (ark->status == TS_STEP_INCOMPLETE) {
103*3a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, X));
104*3a2a065fSHong Zhang       for (j = 0; j < s; j++) w[j] = h * tab->b[j];
105*3a2a065fSHong Zhang       if (ark->is_slow) {
106*3a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_slow, &Xslow));
107*3a2a065fSHong Zhang         PetscCall(VecMAXPY(Xslow, s, w, ark->YdotRHS_slow));
108*3a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_slow, &Xslow));
109*3a2a065fSHong Zhang       }
110*3a2a065fSHong Zhang       if (ark->is_fast) {
111*3a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_fast, &Xfast));
112*3a2a065fSHong Zhang         if (fasthasE) PetscCall(VecMAXPY(Xfast, s, w, ark->YdotRHS_fast));
113*3a2a065fSHong Zhang         for (j = 0; j < s; j++) w[j] = h * tab->bt[j];
114*3a2a065fSHong Zhang         PetscCall(VecMAXPY(Xfast, s, w, ark->YdotI_fast));
115*3a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_fast, &Xfast));
116*3a2a065fSHong Zhang       }
117*3a2a065fSHong Zhang     } else PetscCall(VecCopy(ts->vec_sol, X));
118*3a2a065fSHong Zhang     if (done) *done = PETSC_TRUE;
119*3a2a065fSHong Zhang     PetscFunctionReturn(PETSC_SUCCESS);
120*3a2a065fSHong Zhang   } else if (order == tab->order - 1) {
121*3a2a065fSHong Zhang     if (!tab->bembedt) goto unavailable;
122*3a2a065fSHong Zhang     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
123*3a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, X));
124*3a2a065fSHong Zhang       for (j = 0; j < s; j++) w[j] = h * tab->bembed[j];
125*3a2a065fSHong Zhang       if (ark->is_slow) {
126*3a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_slow, &Xslow));
127*3a2a065fSHong Zhang         PetscCall(VecMAXPY(Xslow, s, w, ark->YdotRHS_slow));
128*3a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_slow, &Xslow));
129*3a2a065fSHong Zhang       }
130*3a2a065fSHong Zhang       if (ark->is_fast) {
131*3a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_fast, &Xfast));
132*3a2a065fSHong Zhang         if (fasthasE) PetscCall(VecMAXPY(Xfast, s, w, ark->YdotRHS_fast));
133*3a2a065fSHong Zhang         for (j = 0; j < s; j++) w[j] = h * tab->bembedt[j];
134*3a2a065fSHong Zhang         PetscCall(VecMAXPY(Xfast, s, w, ark->YdotI_fast));
135*3a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_fast, &Xfast));
136*3a2a065fSHong Zhang       }
137*3a2a065fSHong Zhang     } else { /* Rollback and re-complete using (bet-be,be-b) */
138*3a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, X));
139*3a2a065fSHong Zhang       for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]);
140*3a2a065fSHong Zhang       if (ark->is_slow) {
141*3a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_slow, &Xslow));
142*3a2a065fSHong Zhang         PetscCall(VecMAXPY(Xslow, s, w, ark->YdotRHS_slow));
143*3a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_slow, &Xslow));
144*3a2a065fSHong Zhang       }
145*3a2a065fSHong Zhang       if (ark->is_fast) {
146*3a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_fast, &Xfast));
147*3a2a065fSHong Zhang         if (fasthasE) PetscCall(VecMAXPY(Xfast, tab->s, w, ark->YdotRHS_fast));
148*3a2a065fSHong Zhang         for (j = 0; j < s; j++) w[j] = h * (tab->bembedt[j] - tab->bt[j]);
149*3a2a065fSHong Zhang         PetscCall(VecMAXPY(Xfast, tab->s, w, ark->YdotI_fast));
150*3a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_fast, &Xfast));
151*3a2a065fSHong Zhang       }
152*3a2a065fSHong Zhang     }
153*3a2a065fSHong Zhang     if (done) *done = PETSC_TRUE;
154*3a2a065fSHong Zhang     PetscFunctionReturn(PETSC_SUCCESS);
155*3a2a065fSHong Zhang   }
156*3a2a065fSHong Zhang unavailable:
157*3a2a065fSHong Zhang   if (done) *done = PETSC_FALSE;
158*3a2a065fSHong Zhang   else
159*3a2a065fSHong Zhang     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "ARKIMEX '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT ". Consider using -ts_adapt_type none or a different method that has an embedded estimate.", tab->name,
160*3a2a065fSHong Zhang             tab->order, order);
161*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
162*3a2a065fSHong Zhang }
163*3a2a065fSHong Zhang 
164*3a2a065fSHong Zhang /*
165*3a2a065fSHong Zhang   Additive Runge-Kutta methods for a fast-slow (component-wise partitioned) system in the form
166*3a2a065fSHong Zhang     Ufdot = Ff(t,Uf,Us)
167*3a2a065fSHong Zhang     Usdot = Fs(t,Uf,Us)
168*3a2a065fSHong Zhang 
169*3a2a065fSHong Zhang   Ys[i] = Us_n + dt \sum_{j=1}^{i-1} a[i][j] Fs(t+c_j*dt,Yf[j],Ys[j])
170*3a2a065fSHong Zhang   Ys[i] = Us_n + dt \sum_{j=1}^{i-1} a[i][j] Fs(t+c_j*dt,Yf[j],Ys[j])
171*3a2a065fSHong Zhang 
172*3a2a065fSHong Zhang */
173*3a2a065fSHong Zhang static PetscErrorCode TSStep_ARKIMEX_FastSlowSplit(TS ts)
174*3a2a065fSHong Zhang {
175*3a2a065fSHong Zhang   TS_ARKIMEX      *ark = (TS_ARKIMEX *)ts->data;
176*3a2a065fSHong Zhang   ARKTableau       tab = ark->tableau;
177*3a2a065fSHong Zhang   const PetscInt   s   = tab->s;
178*3a2a065fSHong Zhang   const PetscReal *At = tab->At, *A = tab->A, *ct = tab->ct;
179*3a2a065fSHong Zhang   PetscScalar     *w = ark->work;
180*3a2a065fSHong Zhang   Vec             *Y = ark->Y, Ydot_fast = ark->Ydot, Ydot0_fast = ark->Ydot0, Z = ark->Z, *YdotRHS_fast = ark->YdotRHS_fast, *YdotRHS_slow = ark->YdotRHS_slow, *YdotI_fast = ark->YdotI_fast, Yfast, Yslow, Xfast, Xslow;
181*3a2a065fSHong Zhang   PetscBool        extrapolate = ark->extrapolate;
182*3a2a065fSHong Zhang   TSAdapt          adapt;
183*3a2a065fSHong Zhang   SNES             snes;
184*3a2a065fSHong Zhang   PetscInt         i, j, its, lits;
185*3a2a065fSHong Zhang   PetscInt         rejections = 0;
186*3a2a065fSHong Zhang   PetscBool        fasthasE = PETSC_FALSE, stageok, accept = PETSC_TRUE;
187*3a2a065fSHong Zhang   PetscReal        next_time_step = ts->time_step;
188*3a2a065fSHong Zhang 
189*3a2a065fSHong Zhang   PetscFunctionBegin;
190*3a2a065fSHong Zhang   if (ark->is_fast) PetscCall(TSHasRHSFunction(ark->subts_fast, &fasthasE));
191*3a2a065fSHong Zhang   if (ark->extrapolate && !ark->Y_prev) {
192*3a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_fast, &Xfast));
193*3a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->Y_prev));
194*3a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotI_prev));
195*3a2a065fSHong Zhang     if (fasthasE) PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotRHS_prev));
196*3a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_fast, &Xfast));
197*3a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_slow, &Xslow));
198*3a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_fast, &Xslow));
199*3a2a065fSHong Zhang   }
200*3a2a065fSHong Zhang 
201*3a2a065fSHong Zhang   if (!ts->steprollback) {
202*3a2a065fSHong Zhang     if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */
203*3a2a065fSHong Zhang       PetscCall(VecCopy(YdotI_fast[s - 1], Ydot0_fast));
204*3a2a065fSHong Zhang     }
205*3a2a065fSHong Zhang     if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */
206*3a2a065fSHong Zhang       for (i = 0; i < s; i++) {
207*3a2a065fSHong Zhang         PetscCall(VecISCopy(Y[i], ark->is_fast, SCATTER_REVERSE, ark->Y_prev[i]));
208*3a2a065fSHong Zhang         PetscCall(VecCopy(YdotI_fast[i], ark->YdotI_prev[i]));
209*3a2a065fSHong Zhang         if (fasthasE) PetscCall(VecCopy(YdotRHS_fast[i], ark->YdotRHS_prev[i]));
210*3a2a065fSHong Zhang       }
211*3a2a065fSHong Zhang     }
212*3a2a065fSHong Zhang   }
213*3a2a065fSHong Zhang 
214*3a2a065fSHong Zhang   /* For IMEX we compute a step */
215*3a2a065fSHong Zhang   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) {
216*3a2a065fSHong Zhang     TS ts_start;
217*3a2a065fSHong Zhang     PetscCall(TSClone(ts, &ts_start));
218*3a2a065fSHong Zhang     PetscCall(TSSetSolution(ts_start, ts->vec_sol));
219*3a2a065fSHong Zhang     PetscCall(TSSetTime(ts_start, ts->ptime));
220*3a2a065fSHong Zhang     PetscCall(TSSetMaxSteps(ts_start, ts->steps + 1));
221*3a2a065fSHong Zhang     PetscCall(TSSetMaxTime(ts_start, ts->ptime + ts->time_step));
222*3a2a065fSHong Zhang     PetscCall(TSSetExactFinalTime(ts_start, TS_EXACTFINALTIME_STEPOVER));
223*3a2a065fSHong Zhang     PetscCall(TSSetTimeStep(ts_start, ts->time_step));
224*3a2a065fSHong Zhang     PetscCall(TSSetType(ts_start, TSARKIMEX));
225*3a2a065fSHong Zhang     PetscCall(TSARKIMEXSetFullyImplicit(ts_start, PETSC_TRUE));
226*3a2a065fSHong Zhang     PetscCall(TSARKIMEXSetType(ts_start, TSARKIMEX1BEE));
227*3a2a065fSHong Zhang 
228*3a2a065fSHong Zhang     PetscCall(TSRestartStep(ts_start));
229*3a2a065fSHong Zhang     PetscCall(TSSolve(ts_start, ts->vec_sol));
230*3a2a065fSHong Zhang     PetscCall(TSGetTime(ts_start, &ts->ptime));
231*3a2a065fSHong Zhang     PetscCall(TSGetTimeStep(ts_start, &ts->time_step));
232*3a2a065fSHong Zhang 
233*3a2a065fSHong Zhang     { /* Save the initial slope for the next step */
234*3a2a065fSHong Zhang       TS_ARKIMEX *ark_start = (TS_ARKIMEX *)ts_start->data;
235*3a2a065fSHong Zhang       PetscCall(VecCopy(ark_start->YdotI[ark_start->tableau->s - 1], Ydot0_fast));
236*3a2a065fSHong Zhang     }
237*3a2a065fSHong Zhang     ts->steps++;
238*3a2a065fSHong Zhang     /* Set the correct TS in SNES */
239*3a2a065fSHong Zhang     /* We'll try to bypass this by changing the method on the fly */
240*3a2a065fSHong Zhang     {
241*3a2a065fSHong Zhang       PetscCall(TSGetSNES(ts, &snes));
242*3a2a065fSHong Zhang       PetscCall(TSSetSNES(ts, snes));
243*3a2a065fSHong Zhang     }
244*3a2a065fSHong Zhang     PetscCall(TSDestroy(&ts_start));
245*3a2a065fSHong Zhang   }
246*3a2a065fSHong Zhang 
247*3a2a065fSHong Zhang   ark->status = TS_STEP_INCOMPLETE;
248*3a2a065fSHong Zhang   while (!ts->reason && ark->status != TS_STEP_COMPLETE) {
249*3a2a065fSHong Zhang     PetscReal t = ts->ptime;
250*3a2a065fSHong Zhang     PetscReal h = ts->time_step;
251*3a2a065fSHong Zhang     for (i = 0; i < s; i++) {
252*3a2a065fSHong Zhang       ark->stage_time = t + h * ct[i];
253*3a2a065fSHong Zhang       PetscCall(TSPreStage(ts, ark->stage_time));
254*3a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, Y[i]));
255*3a2a065fSHong Zhang       /* fast components */
256*3a2a065fSHong Zhang       if (ark->is_fast) {
257*3a2a065fSHong Zhang         if (At[i * s + i] == 0) { /* This stage is explicit */
258*3a2a065fSHong Zhang           PetscCheck(i == 0 || ts->equation_type < TS_EQ_IMPLICIT, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Explicit stages other than the first one are not supported for implicit problems");
259*3a2a065fSHong Zhang           PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
260*3a2a065fSHong Zhang           for (j = 0; j < i; j++) w[j] = h * At[i * s + j];
261*3a2a065fSHong Zhang           PetscCall(VecMAXPY(Yfast, i, w, YdotI_fast));
262*3a2a065fSHong Zhang           if (fasthasE) {
263*3a2a065fSHong Zhang             for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
264*3a2a065fSHong Zhang             PetscCall(VecMAXPY(Yfast, i, w, YdotRHS_fast));
265*3a2a065fSHong Zhang           }
266*3a2a065fSHong Zhang           PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
267*3a2a065fSHong Zhang         } else {
268*3a2a065fSHong Zhang           ark->scoeff = 1. / At[i * s + i];
269*3a2a065fSHong Zhang           /* Ydot = shift*(Y-Z) */
270*3a2a065fSHong Zhang           PetscCall(VecISCopy(ts->vec_sol, ark->is_fast, SCATTER_REVERSE, Z));
271*3a2a065fSHong Zhang           for (j = 0; j < i; j++) w[j] = h * At[i * s + j];
272*3a2a065fSHong Zhang           PetscCall(VecMAXPY(Z, i, w, YdotI_fast));
273*3a2a065fSHong Zhang           if (fasthasE) {
274*3a2a065fSHong Zhang             for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
275*3a2a065fSHong Zhang             PetscCall(VecMAXPY(Z, i, w, YdotRHS_fast));
276*3a2a065fSHong Zhang           }
277*3a2a065fSHong Zhang           PetscCall(TSGetSNES(ts, &snes));
278*3a2a065fSHong Zhang           if (ark->is_slow) PetscCall(VecCopy(i > 0 ? Y[i - 1] : ts->vec_sol, ark->Y_snes));
279*3a2a065fSHong Zhang           else ark->Y_snes = Y[i];
280*3a2a065fSHong Zhang           PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
281*3a2a065fSHong Zhang           if (extrapolate && !ts->steprestart) {
282*3a2a065fSHong Zhang             /* Initial guess extrapolated from previous time step stage values */
283*3a2a065fSHong Zhang             PetscCall(TSExtrapolate_ARKIMEX_FastSlowSplit(ts, ct[i], Yfast));
284*3a2a065fSHong Zhang           } else {
285*3a2a065fSHong Zhang             /* Initial guess taken from last stage */
286*3a2a065fSHong Zhang             PetscCall(VecISCopy(i > 0 ? Y[i - 1] : ts->vec_sol, ark->is_fast, SCATTER_REVERSE, Yfast));
287*3a2a065fSHong Zhang           }
288*3a2a065fSHong Zhang           PetscCall(SNESSolve(snes, NULL, Yfast));
289*3a2a065fSHong Zhang           PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
290*3a2a065fSHong Zhang           PetscCall(SNESGetIterationNumber(snes, &its));
291*3a2a065fSHong Zhang           PetscCall(SNESGetLinearSolveIterations(snes, &lits));
292*3a2a065fSHong Zhang           ts->snes_its += its;
293*3a2a065fSHong Zhang           ts->ksp_its += lits;
294*3a2a065fSHong Zhang           PetscCall(TSGetAdapt(ts, &adapt));
295*3a2a065fSHong Zhang           PetscCall(TSAdaptCheckStage(adapt, ts, ark->stage_time, Y[i], &stageok));
296*3a2a065fSHong Zhang           if (!stageok) {
297*3a2a065fSHong Zhang             /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to
298*3a2a065fSHong Zhang              * use extrapolation to initialize the solves on the next attempt. */
299*3a2a065fSHong Zhang             extrapolate = PETSC_FALSE;
300*3a2a065fSHong Zhang             goto reject_step;
301*3a2a065fSHong Zhang           }
302*3a2a065fSHong Zhang         }
303*3a2a065fSHong Zhang 
304*3a2a065fSHong Zhang         if (ts->equation_type >= TS_EQ_IMPLICIT) {
305*3a2a065fSHong Zhang           if (i == 0 && tab->explicit_first_stage) {
306*3a2a065fSHong Zhang             PetscCheck(tab->stiffly_accurate, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "%s %s is not stiffly accurate and therefore explicit-first stage methods cannot be used if the equation is implicit because the slope cannot be evaluated",
307*3a2a065fSHong Zhang                        ((PetscObject)ts)->type_name, ark->tableau->name);
308*3a2a065fSHong Zhang             PetscCall(VecCopy(Ydot0_fast, YdotI_fast[0])); /* YdotI_fast = YdotI_fast(tn-1) */
309*3a2a065fSHong Zhang           } else {
310*3a2a065fSHong Zhang             PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
311*3a2a065fSHong Zhang             PetscCall(VecAXPBYPCZ(YdotI_fast[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Yfast)); /* YdotI = shift*(X-Z) */
312*3a2a065fSHong Zhang             PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
313*3a2a065fSHong Zhang           }
314*3a2a065fSHong Zhang         } else {
315*3a2a065fSHong Zhang           if (i == 0 && tab->explicit_first_stage) {
316*3a2a065fSHong Zhang             PetscCall(VecZeroEntries(Ydot_fast));
317*3a2a065fSHong Zhang             PetscCall(TSComputeIFunction(ark->subts_fast, ark->stage_time, Y[i], Ydot_fast, YdotI_fast[i], ark->imex)); /* YdotI = -G(t,Y,0)   */
318*3a2a065fSHong Zhang             PetscCall(VecScale(YdotI_fast[i], -1.0));
319*3a2a065fSHong Zhang           } else {
320*3a2a065fSHong Zhang             PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
321*3a2a065fSHong Zhang             PetscCall(VecAXPBYPCZ(YdotI_fast[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Yfast)); /* YdotI = shift*(X-Z) */
322*3a2a065fSHong Zhang             PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
323*3a2a065fSHong Zhang           }
324*3a2a065fSHong Zhang           if (fasthasE) PetscCall(TSComputeRHSFunction(ark->subts_fast, ark->stage_time, Y[i], YdotRHS_fast[i]));
325*3a2a065fSHong Zhang         }
326*3a2a065fSHong Zhang       }
327*3a2a065fSHong Zhang       /* slow components */
328*3a2a065fSHong Zhang       if (ark->is_slow) {
329*3a2a065fSHong Zhang         for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
330*3a2a065fSHong Zhang         PetscCall(VecGetSubVector(Y[i], ark->is_slow, &Yslow));
331*3a2a065fSHong Zhang         PetscCall(VecMAXPY(Yslow, i, w, YdotRHS_slow));
332*3a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(Y[i], ark->is_slow, &Yslow));
333*3a2a065fSHong Zhang         PetscCall(TSComputeRHSFunction(ark->subts_slow, ark->stage_time, Y[i], YdotRHS_slow[i]));
334*3a2a065fSHong Zhang       }
335*3a2a065fSHong Zhang       PetscCall(TSPostStage(ts, ark->stage_time, i, Y));
336*3a2a065fSHong Zhang     }
337*3a2a065fSHong Zhang     ark->status = TS_STEP_INCOMPLETE;
338*3a2a065fSHong Zhang     PetscCall(TSEvaluateStep_ARKIMEX_FastSlowSplit(ts, tab->order, ts->vec_sol, NULL));
339*3a2a065fSHong Zhang     ark->status = TS_STEP_PENDING;
340*3a2a065fSHong Zhang     PetscCall(TSGetAdapt(ts, &adapt));
341*3a2a065fSHong Zhang     PetscCall(TSAdaptCandidatesClear(adapt));
342*3a2a065fSHong Zhang     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
343*3a2a065fSHong Zhang     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
344*3a2a065fSHong Zhang     ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
345*3a2a065fSHong Zhang     if (!accept) { /* Roll back the current step */
346*3a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol0, ts->vec_sol));
347*3a2a065fSHong Zhang       ts->time_step = next_time_step;
348*3a2a065fSHong Zhang       goto reject_step;
349*3a2a065fSHong Zhang     }
350*3a2a065fSHong Zhang 
351*3a2a065fSHong Zhang     ts->ptime += ts->time_step;
352*3a2a065fSHong Zhang     ts->time_step = next_time_step;
353*3a2a065fSHong Zhang     break;
354*3a2a065fSHong Zhang 
355*3a2a065fSHong Zhang   reject_step:
356*3a2a065fSHong Zhang     ts->reject++;
357*3a2a065fSHong Zhang     accept = PETSC_FALSE;
358*3a2a065fSHong Zhang     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
359*3a2a065fSHong Zhang       ts->reason = TS_DIVERGED_STEP_REJECTED;
360*3a2a065fSHong Zhang       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
361*3a2a065fSHong Zhang     }
362*3a2a065fSHong Zhang   }
363*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
364*3a2a065fSHong Zhang }
365*3a2a065fSHong Zhang 
366*3a2a065fSHong Zhang static PetscErrorCode TSSetUp_ARKIMEX_FastSlowSplit(TS ts)
367*3a2a065fSHong Zhang {
368*3a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
369*3a2a065fSHong Zhang   ARKTableau  tab = ark->tableau;
370*3a2a065fSHong Zhang   Vec         Xfast, Xslow;
371*3a2a065fSHong Zhang 
372*3a2a065fSHong Zhang   PetscFunctionBegin;
373*3a2a065fSHong Zhang   PetscCall(PetscMalloc1(2 * tab->s, &ark->work));
374*3a2a065fSHong Zhang   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y));
375*3a2a065fSHong Zhang   PetscCall(TSRHSSplitGetIS(ts, "slow", &ark->is_slow));
376*3a2a065fSHong Zhang   PetscCall(TSRHSSplitGetIS(ts, "fast", &ark->is_fast));
377*3a2a065fSHong Zhang   PetscCheck(ark->is_slow || ark->is_fast, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' or 'fast' or both in order to use -ts_arkimex_fastslow true");
378*3a2a065fSHong Zhang   /* The following vectors need to be resized */
379*3a2a065fSHong Zhang   PetscCall(VecDestroy(&ark->Ydot));
380*3a2a065fSHong Zhang   PetscCall(VecDestroy(&ark->Ydot0));
381*3a2a065fSHong Zhang   PetscCall(VecDestroy(&ark->Z));
382*3a2a065fSHong Zhang   PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_fast));
383*3a2a065fSHong Zhang   if (ark->extrapolate && ark->is_slow) { // need to resize these vectors if the fast subvectors is smaller than their original counterparts (which means IS)
384*3a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->Y_prev));
385*3a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_prev));
386*3a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_prev));
387*3a2a065fSHong Zhang   }
388*3a2a065fSHong Zhang   /* Allocate working vectors */
389*3a2a065fSHong Zhang   if (ark->is_fast && ark->is_slow) PetscCall(VecDuplicate(ts->vec_sol, &ark->Y_snes)); // need an additional work vector for SNES
390*3a2a065fSHong Zhang   if (ark->is_fast) {
391*3a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_fast, &Xfast));
392*3a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotRHS_fast));
393*3a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotI_fast));
394*3a2a065fSHong Zhang     PetscCall(VecDuplicate(Xfast, &ark->Ydot));
395*3a2a065fSHong Zhang     PetscCall(VecDuplicate(Xfast, &ark->Ydot0));
396*3a2a065fSHong Zhang     PetscCall(VecDuplicate(Xfast, &ark->Z));
397*3a2a065fSHong Zhang     if (ark->extrapolate) {
398*3a2a065fSHong Zhang       PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->Y_prev));
399*3a2a065fSHong Zhang       PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotI_prev));
400*3a2a065fSHong Zhang       PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotRHS_prev));
401*3a2a065fSHong Zhang     }
402*3a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_fast, &Xfast));
403*3a2a065fSHong Zhang   }
404*3a2a065fSHong Zhang   if (ark->is_slow) {
405*3a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_slow, &Xslow));
406*3a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xslow, tab->s, &ark->YdotRHS_slow));
407*3a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_slow, &Xslow));
408*3a2a065fSHong Zhang   }
409*3a2a065fSHong Zhang   ts->ops->step         = TSStep_ARKIMEX_FastSlowSplit;
410*3a2a065fSHong Zhang   ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX_FastSlowSplit;
411*3a2a065fSHong Zhang   PetscCall(TSARKIMEXSetSplits(ts));
412*3a2a065fSHong Zhang   if (ark->subts_fast) { // subts SNESJacobian is set when users set the subts Jacobian, but the main ts SNESJacobian needs to be set too
413*3a2a065fSHong Zhang     SNES snes, snes_fast;
414*3a2a065fSHong Zhang     PetscErrorCode (*func)(SNES, Vec, Mat, Mat, void *);
415*3a2a065fSHong Zhang     ts->matchsnesdm = PETSC_FALSE;
416*3a2a065fSHong Zhang     PetscCall(TSGetSNES(ts, &snes));
417*3a2a065fSHong Zhang     PetscCall(TSGetSNES(ark->subts_fast, &snes_fast));
418*3a2a065fSHong Zhang     PetscCall(SNESGetJacobian(snes_fast, NULL, NULL, &func, NULL));
419*3a2a065fSHong Zhang     if (func == SNESTSFormJacobian) PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESTSFormJacobian, ts));
420*3a2a065fSHong Zhang   }
421*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
422*3a2a065fSHong Zhang }
423*3a2a065fSHong Zhang 
424*3a2a065fSHong Zhang static PetscErrorCode TSReset_ARKIMEX_FastSlowSplit(TS ts)
425*3a2a065fSHong Zhang {
426*3a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
427*3a2a065fSHong Zhang   ARKTableau  tab = ark->tableau;
428*3a2a065fSHong Zhang 
429*3a2a065fSHong Zhang   PetscFunctionBegin;
430*3a2a065fSHong Zhang   if (tab) {
431*3a2a065fSHong Zhang     PetscCall(PetscFree(ark->work));
432*3a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->Y));
433*3a2a065fSHong Zhang     if (ark->is_fast && ark->is_slow) PetscCall(VecDestroy(&ark->Y_snes));
434*3a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_slow));
435*3a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_fast));
436*3a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_fast));
437*3a2a065fSHong Zhang     PetscCall(VecDestroy(&ark->Ydot));
438*3a2a065fSHong Zhang     PetscCall(VecDestroy(&ark->Ydot0));
439*3a2a065fSHong Zhang     PetscCall(VecDestroy(&ark->Z));
440*3a2a065fSHong Zhang     if (ark->extrapolate) {
441*3a2a065fSHong Zhang       PetscCall(VecDestroyVecs(tab->s, &ark->Y_prev));
442*3a2a065fSHong Zhang       PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_prev));
443*3a2a065fSHong Zhang       PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_prev));
444*3a2a065fSHong Zhang     }
445*3a2a065fSHong Zhang   }
446*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
447*3a2a065fSHong Zhang }
448*3a2a065fSHong Zhang 
449*3a2a065fSHong Zhang PetscErrorCode TSARKIMEXSetFastSlowSplit_ARKIMEX(TS ts, PetscBool fastslowsplit)
450*3a2a065fSHong Zhang {
451*3a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
452*3a2a065fSHong Zhang 
453*3a2a065fSHong Zhang   PetscFunctionBegin;
454*3a2a065fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
455*3a2a065fSHong Zhang   ark->fastslowsplit = fastslowsplit;
456*3a2a065fSHong Zhang   if (fastslowsplit) {
457*3a2a065fSHong Zhang     ts->matchsnesdm = PETSC_FALSE;
458*3a2a065fSHong Zhang     PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_ARKIMEX_FastSlowSplit_C", TSSetUp_ARKIMEX_FastSlowSplit));
459*3a2a065fSHong Zhang     PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_ARKIMEX_FastSlowSplit_C", TSReset_ARKIMEX_FastSlowSplit));
460*3a2a065fSHong Zhang   }
461*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
462*3a2a065fSHong Zhang }
463*3a2a065fSHong Zhang 
464*3a2a065fSHong Zhang PetscErrorCode TSARKIMEXGetFastSlowSplit_ARKIMEX(TS ts, PetscBool *fastslowsplit)
465*3a2a065fSHong Zhang {
466*3a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
467*3a2a065fSHong Zhang 
468*3a2a065fSHong Zhang   PetscFunctionBegin;
469*3a2a065fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
470*3a2a065fSHong Zhang   *fastslowsplit = ark->fastslowsplit;
471*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
472*3a2a065fSHong Zhang }
473