xref: /petsc/src/ts/impls/arkimex/fsarkimex.c (revision 4bd3aaa3b0a0aba9d82762bb2193c2700c26f8e1)
13a2a065fSHong Zhang #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
23a2a065fSHong Zhang #include <petscdm.h>
33a2a065fSHong Zhang #include <../src/ts/impls/arkimex/arkimex.h>
43a2a065fSHong Zhang #include <../src/ts/impls/arkimex/fsarkimex.h>
53a2a065fSHong Zhang 
63a2a065fSHong Zhang static PetscErrorCode TSARKIMEXSetSplits(TS ts)
73a2a065fSHong Zhang {
83a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
93a2a065fSHong Zhang   DM          dm, subdm, newdm;
103a2a065fSHong Zhang 
113a2a065fSHong Zhang   PetscFunctionBegin;
123a2a065fSHong Zhang   PetscCall(TSRHSSplitGetSubTS(ts, "slow", &ark->subts_slow));
133a2a065fSHong Zhang   PetscCall(TSRHSSplitGetSubTS(ts, "fast", &ark->subts_fast));
143a2a065fSHong Zhang   /* Only copy the DM */
153a2a065fSHong Zhang   PetscCall(TSGetDM(ts, &dm));
163a2a065fSHong Zhang   if (ark->subts_slow) {
173a2a065fSHong Zhang     PetscCall(DMClone(dm, &newdm));
183a2a065fSHong Zhang     PetscCall(TSGetDM(ark->subts_slow, &subdm));
193a2a065fSHong Zhang     PetscCall(DMCopyDMTS(subdm, newdm));
203a2a065fSHong Zhang     PetscCall(TSSetDM(ark->subts_slow, newdm));
213a2a065fSHong Zhang     PetscCall(DMDestroy(&newdm));
223a2a065fSHong Zhang   }
233a2a065fSHong Zhang   if (ark->subts_fast) {
243a2a065fSHong Zhang     PetscCall(DMClone(dm, &newdm));
253a2a065fSHong Zhang     PetscCall(TSGetDM(ark->subts_fast, &subdm));
263a2a065fSHong Zhang     PetscCall(DMCopyDMTS(subdm, newdm));
273a2a065fSHong Zhang     PetscCall(TSSetDM(ark->subts_fast, newdm));
283a2a065fSHong Zhang     PetscCall(DMDestroy(&newdm));
293a2a065fSHong Zhang   }
303a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
313a2a065fSHong Zhang }
323a2a065fSHong Zhang 
333a2a065fSHong Zhang static PetscErrorCode TSExtrapolate_ARKIMEX_FastSlowSplit(TS ts, PetscReal c, Vec X)
343a2a065fSHong Zhang {
353a2a065fSHong Zhang   TS_ARKIMEX      *ark = (TS_ARKIMEX *)ts->data;
363a2a065fSHong Zhang   ARKTableau       tab = ark->tableau;
373a2a065fSHong Zhang   PetscInt         s = tab->s, pinterp = tab->pinterp, i, j;
383a2a065fSHong Zhang   PetscReal        h, h_prev, t, tt;
393a2a065fSHong Zhang   PetscScalar     *bt = ark->work, *b = ark->work + s;
403a2a065fSHong Zhang   const PetscReal *Bt = tab->binterpt, *B = tab->binterp;
413a2a065fSHong Zhang   PetscBool        fasthasE;
423a2a065fSHong Zhang 
433a2a065fSHong Zhang   PetscFunctionBegin;
443a2a065fSHong Zhang   PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name);
453a2a065fSHong Zhang   h      = ts->time_step;
463a2a065fSHong Zhang   h_prev = ts->ptime - ts->ptime_prev;
473a2a065fSHong Zhang   t      = 1 + h / h_prev * c;
483a2a065fSHong Zhang   for (i = 0; i < s; i++) bt[i] = b[i] = 0;
493a2a065fSHong Zhang   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
503a2a065fSHong Zhang     for (i = 0; i < s; i++) {
513a2a065fSHong Zhang       bt[i] += h * Bt[i * pinterp + j] * tt;
523a2a065fSHong Zhang       b[i] += h * B[i * pinterp + j] * tt;
533a2a065fSHong Zhang     }
543a2a065fSHong Zhang   }
553a2a065fSHong Zhang   PetscCheck(ark->Y_prev, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Stages from previous step have not been stored");
563a2a065fSHong Zhang   PetscCall(VecCopy(ark->Y_prev[0], X));
573a2a065fSHong Zhang   PetscCall(VecMAXPY(X, s, bt, ark->YdotI_prev));
583a2a065fSHong Zhang   PetscCall(TSHasRHSFunction(ark->subts_fast, &fasthasE));
593a2a065fSHong Zhang   if (fasthasE) PetscCall(VecMAXPY(X, s, b, ark->YdotRHS_prev));
603a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
613a2a065fSHong Zhang }
623a2a065fSHong Zhang 
633a2a065fSHong Zhang /*
643a2a065fSHong Zhang  The step completion formula is
653a2a065fSHong Zhang 
663a2a065fSHong Zhang  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
673a2a065fSHong Zhang 
683a2a065fSHong Zhang  This function can be called before or after ts->vec_sol has been updated.
693a2a065fSHong Zhang  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
703a2a065fSHong Zhang  We can write
713a2a065fSHong Zhang 
723a2a065fSHong Zhang  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
733a2a065fSHong Zhang      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
743a2a065fSHong Zhang      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
753a2a065fSHong Zhang 
763a2a065fSHong Zhang  so we can evaluate the method with different order even after the step has been optimistically completed.
773a2a065fSHong Zhang */
783a2a065fSHong Zhang static PetscErrorCode TSEvaluateStep_ARKIMEX_FastSlowSplit(TS ts, PetscInt order, Vec X, PetscBool *done)
793a2a065fSHong Zhang {
803a2a065fSHong Zhang   TS_ARKIMEX  *ark = (TS_ARKIMEX *)ts->data;
813a2a065fSHong Zhang   ARKTableau   tab = ark->tableau;
823a2a065fSHong Zhang   Vec          Xfast, Xslow;
833a2a065fSHong Zhang   PetscScalar *w = ark->work;
843a2a065fSHong Zhang   PetscReal    h;
853a2a065fSHong Zhang   PetscInt     s = tab->s, j;
863a2a065fSHong Zhang   PetscBool    fasthasE;
873a2a065fSHong Zhang 
883a2a065fSHong Zhang   PetscFunctionBegin;
893a2a065fSHong Zhang   switch (ark->status) {
903a2a065fSHong Zhang   case TS_STEP_INCOMPLETE:
913a2a065fSHong Zhang   case TS_STEP_PENDING:
923a2a065fSHong Zhang     h = ts->time_step;
933a2a065fSHong Zhang     break;
943a2a065fSHong Zhang   case TS_STEP_COMPLETE:
953a2a065fSHong Zhang     h = ts->ptime - ts->ptime_prev;
963a2a065fSHong Zhang     break;
973a2a065fSHong Zhang   default:
983a2a065fSHong Zhang     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
993a2a065fSHong Zhang   }
1003a2a065fSHong Zhang   if (ark->is_fast) PetscCall(TSHasRHSFunction(ark->subts_fast, &fasthasE));
1013a2a065fSHong Zhang   if (order == tab->order) {
1023a2a065fSHong Zhang     if (ark->status == TS_STEP_INCOMPLETE) {
1033a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, X));
1043a2a065fSHong Zhang       for (j = 0; j < s; j++) w[j] = h * tab->b[j];
1053a2a065fSHong Zhang       if (ark->is_slow) {
1063a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_slow, &Xslow));
1073a2a065fSHong Zhang         PetscCall(VecMAXPY(Xslow, s, w, ark->YdotRHS_slow));
1083a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_slow, &Xslow));
1093a2a065fSHong Zhang       }
1103a2a065fSHong Zhang       if (ark->is_fast) {
1113a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_fast, &Xfast));
1123a2a065fSHong Zhang         if (fasthasE) PetscCall(VecMAXPY(Xfast, s, w, ark->YdotRHS_fast));
1133a2a065fSHong Zhang         for (j = 0; j < s; j++) w[j] = h * tab->bt[j];
1143a2a065fSHong Zhang         PetscCall(VecMAXPY(Xfast, s, w, ark->YdotI_fast));
1153a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_fast, &Xfast));
1163a2a065fSHong Zhang       }
1173a2a065fSHong Zhang     } else PetscCall(VecCopy(ts->vec_sol, X));
1183a2a065fSHong Zhang     if (done) *done = PETSC_TRUE;
1193a2a065fSHong Zhang     PetscFunctionReturn(PETSC_SUCCESS);
1203a2a065fSHong Zhang   } else if (order == tab->order - 1) {
1213a2a065fSHong Zhang     if (!tab->bembedt) goto unavailable;
1223a2a065fSHong Zhang     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
1233a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, X));
1243a2a065fSHong Zhang       for (j = 0; j < s; j++) w[j] = h * tab->bembed[j];
1253a2a065fSHong Zhang       if (ark->is_slow) {
1263a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_slow, &Xslow));
1273a2a065fSHong Zhang         PetscCall(VecMAXPY(Xslow, s, w, ark->YdotRHS_slow));
1283a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_slow, &Xslow));
1293a2a065fSHong Zhang       }
1303a2a065fSHong Zhang       if (ark->is_fast) {
1313a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_fast, &Xfast));
1323a2a065fSHong Zhang         if (fasthasE) PetscCall(VecMAXPY(Xfast, s, w, ark->YdotRHS_fast));
1333a2a065fSHong Zhang         for (j = 0; j < s; j++) w[j] = h * tab->bembedt[j];
1343a2a065fSHong Zhang         PetscCall(VecMAXPY(Xfast, s, w, ark->YdotI_fast));
1353a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_fast, &Xfast));
1363a2a065fSHong Zhang       }
1373a2a065fSHong Zhang     } else { /* Rollback and re-complete using (bet-be,be-b) */
1383a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, X));
1393a2a065fSHong Zhang       for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]);
1403a2a065fSHong Zhang       if (ark->is_slow) {
1413a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_slow, &Xslow));
1423a2a065fSHong Zhang         PetscCall(VecMAXPY(Xslow, s, w, ark->YdotRHS_slow));
1433a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_slow, &Xslow));
1443a2a065fSHong Zhang       }
1453a2a065fSHong Zhang       if (ark->is_fast) {
1463a2a065fSHong Zhang         PetscCall(VecGetSubVector(X, ark->is_fast, &Xfast));
1473a2a065fSHong Zhang         if (fasthasE) PetscCall(VecMAXPY(Xfast, tab->s, w, ark->YdotRHS_fast));
1483a2a065fSHong Zhang         for (j = 0; j < s; j++) w[j] = h * (tab->bembedt[j] - tab->bt[j]);
1493a2a065fSHong Zhang         PetscCall(VecMAXPY(Xfast, tab->s, w, ark->YdotI_fast));
1503a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(X, ark->is_fast, &Xfast));
1513a2a065fSHong Zhang       }
1523a2a065fSHong Zhang     }
1533a2a065fSHong Zhang     if (done) *done = PETSC_TRUE;
1543a2a065fSHong Zhang     PetscFunctionReturn(PETSC_SUCCESS);
1553a2a065fSHong Zhang   }
1563a2a065fSHong Zhang unavailable:
1573a2a065fSHong Zhang   if (done) *done = PETSC_FALSE;
1583a2a065fSHong Zhang   else
1593a2a065fSHong 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,
1603a2a065fSHong Zhang             tab->order, order);
1613a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
1623a2a065fSHong Zhang }
1633a2a065fSHong Zhang 
1643a2a065fSHong Zhang /*
1653a2a065fSHong Zhang   Additive Runge-Kutta methods for a fast-slow (component-wise partitioned) system in the form
1663a2a065fSHong Zhang     Ufdot = Ff(t,Uf,Us)
1673a2a065fSHong Zhang     Usdot = Fs(t,Uf,Us)
1683a2a065fSHong Zhang 
1693a2a065fSHong Zhang   Ys[i] = Us_n + dt \sum_{j=1}^{i-1} a[i][j] Fs(t+c_j*dt,Yf[j],Ys[j])
1703a2a065fSHong Zhang   Ys[i] = Us_n + dt \sum_{j=1}^{i-1} a[i][j] Fs(t+c_j*dt,Yf[j],Ys[j])
1713a2a065fSHong Zhang 
1723a2a065fSHong Zhang */
1733a2a065fSHong Zhang static PetscErrorCode TSStep_ARKIMEX_FastSlowSplit(TS ts)
1743a2a065fSHong Zhang {
1753a2a065fSHong Zhang   TS_ARKIMEX      *ark = (TS_ARKIMEX *)ts->data;
1763a2a065fSHong Zhang   ARKTableau       tab = ark->tableau;
1773a2a065fSHong Zhang   const PetscInt   s   = tab->s;
1783a2a065fSHong Zhang   const PetscReal *At = tab->At, *A = tab->A, *ct = tab->ct;
1793a2a065fSHong Zhang   PetscScalar     *w = ark->work;
1803a2a065fSHong 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;
1813a2a065fSHong Zhang   PetscBool        extrapolate = ark->extrapolate;
1823a2a065fSHong Zhang   TSAdapt          adapt;
1833a2a065fSHong Zhang   SNES             snes;
1843a2a065fSHong Zhang   PetscInt         i, j, its, lits;
1853a2a065fSHong Zhang   PetscInt         rejections = 0;
1863a2a065fSHong Zhang   PetscBool        fasthasE = PETSC_FALSE, stageok, accept = PETSC_TRUE;
1873a2a065fSHong Zhang   PetscReal        next_time_step = ts->time_step;
1883a2a065fSHong Zhang 
1893a2a065fSHong Zhang   PetscFunctionBegin;
1903a2a065fSHong Zhang   if (ark->is_fast) PetscCall(TSHasRHSFunction(ark->subts_fast, &fasthasE));
1913a2a065fSHong Zhang   if (ark->extrapolate && !ark->Y_prev) {
1923a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_fast, &Xfast));
1933a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->Y_prev));
1943a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotI_prev));
1953a2a065fSHong Zhang     if (fasthasE) PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotRHS_prev));
1963a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_fast, &Xfast));
1973a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_slow, &Xslow));
1983a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_fast, &Xslow));
1993a2a065fSHong Zhang   }
2003a2a065fSHong Zhang 
2013a2a065fSHong Zhang   if (!ts->steprollback) {
2023a2a065fSHong Zhang     if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */
2033a2a065fSHong Zhang       PetscCall(VecCopy(YdotI_fast[s - 1], Ydot0_fast));
2043a2a065fSHong Zhang     }
2053a2a065fSHong Zhang     if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */
2063a2a065fSHong Zhang       for (i = 0; i < s; i++) {
2073a2a065fSHong Zhang         PetscCall(VecISCopy(Y[i], ark->is_fast, SCATTER_REVERSE, ark->Y_prev[i]));
2083a2a065fSHong Zhang         PetscCall(VecCopy(YdotI_fast[i], ark->YdotI_prev[i]));
2093a2a065fSHong Zhang         if (fasthasE) PetscCall(VecCopy(YdotRHS_fast[i], ark->YdotRHS_prev[i]));
2103a2a065fSHong Zhang       }
2113a2a065fSHong Zhang     }
2123a2a065fSHong Zhang   }
2133a2a065fSHong Zhang 
2143a2a065fSHong Zhang   /* For IMEX we compute a step */
2153a2a065fSHong Zhang   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) {
2163a2a065fSHong Zhang     TS ts_start;
2173a2a065fSHong Zhang     PetscCall(TSClone(ts, &ts_start));
2183a2a065fSHong Zhang     PetscCall(TSSetSolution(ts_start, ts->vec_sol));
2193a2a065fSHong Zhang     PetscCall(TSSetTime(ts_start, ts->ptime));
2203a2a065fSHong Zhang     PetscCall(TSSetMaxSteps(ts_start, ts->steps + 1));
2213a2a065fSHong Zhang     PetscCall(TSSetMaxTime(ts_start, ts->ptime + ts->time_step));
2223a2a065fSHong Zhang     PetscCall(TSSetExactFinalTime(ts_start, TS_EXACTFINALTIME_STEPOVER));
2233a2a065fSHong Zhang     PetscCall(TSSetTimeStep(ts_start, ts->time_step));
2243a2a065fSHong Zhang     PetscCall(TSSetType(ts_start, TSARKIMEX));
2253a2a065fSHong Zhang     PetscCall(TSARKIMEXSetFullyImplicit(ts_start, PETSC_TRUE));
2263a2a065fSHong Zhang     PetscCall(TSARKIMEXSetType(ts_start, TSARKIMEX1BEE));
2273a2a065fSHong Zhang 
2283a2a065fSHong Zhang     PetscCall(TSRestartStep(ts_start));
2293a2a065fSHong Zhang     PetscCall(TSSolve(ts_start, ts->vec_sol));
2303a2a065fSHong Zhang     PetscCall(TSGetTime(ts_start, &ts->ptime));
2313a2a065fSHong Zhang     PetscCall(TSGetTimeStep(ts_start, &ts->time_step));
2323a2a065fSHong Zhang 
2333a2a065fSHong Zhang     { /* Save the initial slope for the next step */
2343a2a065fSHong Zhang       TS_ARKIMEX *ark_start = (TS_ARKIMEX *)ts_start->data;
2353a2a065fSHong Zhang       PetscCall(VecCopy(ark_start->YdotI[ark_start->tableau->s - 1], Ydot0_fast));
2363a2a065fSHong Zhang     }
2373a2a065fSHong Zhang     ts->steps++;
2383a2a065fSHong Zhang     /* Set the correct TS in SNES */
2393a2a065fSHong Zhang     /* We'll try to bypass this by changing the method on the fly */
2403a2a065fSHong Zhang     {
241*4bd3aaa3SHong Zhang       PetscCall(TSRHSSplitGetSNES(ts, &snes));
242*4bd3aaa3SHong Zhang       PetscCall(TSRHSSplitSetSNES(ts, snes));
2433a2a065fSHong Zhang     }
2443a2a065fSHong Zhang     PetscCall(TSDestroy(&ts_start));
2453a2a065fSHong Zhang   }
2463a2a065fSHong Zhang 
2473a2a065fSHong Zhang   ark->status = TS_STEP_INCOMPLETE;
2483a2a065fSHong Zhang   while (!ts->reason && ark->status != TS_STEP_COMPLETE) {
2493a2a065fSHong Zhang     PetscReal t = ts->ptime;
2503a2a065fSHong Zhang     PetscReal h = ts->time_step;
2513a2a065fSHong Zhang     for (i = 0; i < s; i++) {
2523a2a065fSHong Zhang       ark->stage_time = t + h * ct[i];
2533a2a065fSHong Zhang       PetscCall(TSPreStage(ts, ark->stage_time));
2543a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol, Y[i]));
2553a2a065fSHong Zhang       /* fast components */
2563a2a065fSHong Zhang       if (ark->is_fast) {
2573a2a065fSHong Zhang         if (At[i * s + i] == 0) { /* This stage is explicit */
2583a2a065fSHong 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");
2593a2a065fSHong Zhang           PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
2603a2a065fSHong Zhang           for (j = 0; j < i; j++) w[j] = h * At[i * s + j];
2613a2a065fSHong Zhang           PetscCall(VecMAXPY(Yfast, i, w, YdotI_fast));
2623a2a065fSHong Zhang           if (fasthasE) {
2633a2a065fSHong Zhang             for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
2643a2a065fSHong Zhang             PetscCall(VecMAXPY(Yfast, i, w, YdotRHS_fast));
2653a2a065fSHong Zhang           }
2663a2a065fSHong Zhang           PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
2673a2a065fSHong Zhang         } else {
2683a2a065fSHong Zhang           ark->scoeff = 1. / At[i * s + i];
2693a2a065fSHong Zhang           /* Ydot = shift*(Y-Z) */
2703a2a065fSHong Zhang           PetscCall(VecISCopy(ts->vec_sol, ark->is_fast, SCATTER_REVERSE, Z));
2713a2a065fSHong Zhang           for (j = 0; j < i; j++) w[j] = h * At[i * s + j];
2723a2a065fSHong Zhang           PetscCall(VecMAXPY(Z, i, w, YdotI_fast));
2733a2a065fSHong Zhang           if (fasthasE) {
2743a2a065fSHong Zhang             for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
2753a2a065fSHong Zhang             PetscCall(VecMAXPY(Z, i, w, YdotRHS_fast));
2763a2a065fSHong Zhang           }
277*4bd3aaa3SHong Zhang           PetscCall(TSRHSSplitGetSNES(ts, &snes));
2783a2a065fSHong Zhang           if (ark->is_slow) PetscCall(VecCopy(i > 0 ? Y[i - 1] : ts->vec_sol, ark->Y_snes));
2793a2a065fSHong Zhang           else ark->Y_snes = Y[i];
2803a2a065fSHong Zhang           PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
2813a2a065fSHong Zhang           if (extrapolate && !ts->steprestart) {
2823a2a065fSHong Zhang             /* Initial guess extrapolated from previous time step stage values */
2833a2a065fSHong Zhang             PetscCall(TSExtrapolate_ARKIMEX_FastSlowSplit(ts, ct[i], Yfast));
2843a2a065fSHong Zhang           } else {
2853a2a065fSHong Zhang             /* Initial guess taken from last stage */
2863a2a065fSHong Zhang             PetscCall(VecISCopy(i > 0 ? Y[i - 1] : ts->vec_sol, ark->is_fast, SCATTER_REVERSE, Yfast));
2873a2a065fSHong Zhang           }
2883a2a065fSHong Zhang           PetscCall(SNESSolve(snes, NULL, Yfast));
2893a2a065fSHong Zhang           PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
2903a2a065fSHong Zhang           PetscCall(SNESGetIterationNumber(snes, &its));
2913a2a065fSHong Zhang           PetscCall(SNESGetLinearSolveIterations(snes, &lits));
2923a2a065fSHong Zhang           ts->snes_its += its;
2933a2a065fSHong Zhang           ts->ksp_its += lits;
2943a2a065fSHong Zhang           PetscCall(TSGetAdapt(ts, &adapt));
2953a2a065fSHong Zhang           PetscCall(TSAdaptCheckStage(adapt, ts, ark->stage_time, Y[i], &stageok));
2963a2a065fSHong Zhang           if (!stageok) {
2973a2a065fSHong Zhang             /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to
2983a2a065fSHong Zhang              * use extrapolation to initialize the solves on the next attempt. */
2993a2a065fSHong Zhang             extrapolate = PETSC_FALSE;
3003a2a065fSHong Zhang             goto reject_step;
3013a2a065fSHong Zhang           }
3023a2a065fSHong Zhang         }
3033a2a065fSHong Zhang 
3043a2a065fSHong Zhang         if (ts->equation_type >= TS_EQ_IMPLICIT) {
3053a2a065fSHong Zhang           if (i == 0 && tab->explicit_first_stage) {
3063a2a065fSHong 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",
3073a2a065fSHong Zhang                        ((PetscObject)ts)->type_name, ark->tableau->name);
3083a2a065fSHong Zhang             PetscCall(VecCopy(Ydot0_fast, YdotI_fast[0])); /* YdotI_fast = YdotI_fast(tn-1) */
3093a2a065fSHong Zhang           } else {
3103a2a065fSHong Zhang             PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
3113a2a065fSHong Zhang             PetscCall(VecAXPBYPCZ(YdotI_fast[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Yfast)); /* YdotI = shift*(X-Z) */
3123a2a065fSHong Zhang             PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
3133a2a065fSHong Zhang           }
3143a2a065fSHong Zhang         } else {
3153a2a065fSHong Zhang           if (i == 0 && tab->explicit_first_stage) {
3163a2a065fSHong Zhang             PetscCall(VecZeroEntries(Ydot_fast));
3173a2a065fSHong Zhang             PetscCall(TSComputeIFunction(ark->subts_fast, ark->stage_time, Y[i], Ydot_fast, YdotI_fast[i], ark->imex)); /* YdotI = -G(t,Y,0)   */
3183a2a065fSHong Zhang             PetscCall(VecScale(YdotI_fast[i], -1.0));
3193a2a065fSHong Zhang           } else {
3203a2a065fSHong Zhang             PetscCall(VecGetSubVector(Y[i], ark->is_fast, &Yfast));
3213a2a065fSHong Zhang             PetscCall(VecAXPBYPCZ(YdotI_fast[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Yfast)); /* YdotI = shift*(X-Z) */
3223a2a065fSHong Zhang             PetscCall(VecRestoreSubVector(Y[i], ark->is_fast, &Yfast));
3233a2a065fSHong Zhang           }
3243a2a065fSHong Zhang           if (fasthasE) PetscCall(TSComputeRHSFunction(ark->subts_fast, ark->stage_time, Y[i], YdotRHS_fast[i]));
3253a2a065fSHong Zhang         }
3263a2a065fSHong Zhang       }
3273a2a065fSHong Zhang       /* slow components */
3283a2a065fSHong Zhang       if (ark->is_slow) {
3293a2a065fSHong Zhang         for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
3303a2a065fSHong Zhang         PetscCall(VecGetSubVector(Y[i], ark->is_slow, &Yslow));
3313a2a065fSHong Zhang         PetscCall(VecMAXPY(Yslow, i, w, YdotRHS_slow));
3323a2a065fSHong Zhang         PetscCall(VecRestoreSubVector(Y[i], ark->is_slow, &Yslow));
3333a2a065fSHong Zhang         PetscCall(TSComputeRHSFunction(ark->subts_slow, ark->stage_time, Y[i], YdotRHS_slow[i]));
3343a2a065fSHong Zhang       }
3353a2a065fSHong Zhang       PetscCall(TSPostStage(ts, ark->stage_time, i, Y));
3363a2a065fSHong Zhang     }
3373a2a065fSHong Zhang     ark->status = TS_STEP_INCOMPLETE;
3383a2a065fSHong Zhang     PetscCall(TSEvaluateStep_ARKIMEX_FastSlowSplit(ts, tab->order, ts->vec_sol, NULL));
3393a2a065fSHong Zhang     ark->status = TS_STEP_PENDING;
3403a2a065fSHong Zhang     PetscCall(TSGetAdapt(ts, &adapt));
3413a2a065fSHong Zhang     PetscCall(TSAdaptCandidatesClear(adapt));
3423a2a065fSHong Zhang     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
3433a2a065fSHong Zhang     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
3443a2a065fSHong Zhang     ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
3453a2a065fSHong Zhang     if (!accept) { /* Roll back the current step */
3463a2a065fSHong Zhang       PetscCall(VecCopy(ts->vec_sol0, ts->vec_sol));
3473a2a065fSHong Zhang       ts->time_step = next_time_step;
3483a2a065fSHong Zhang       goto reject_step;
3493a2a065fSHong Zhang     }
3503a2a065fSHong Zhang 
3513a2a065fSHong Zhang     ts->ptime += ts->time_step;
3523a2a065fSHong Zhang     ts->time_step = next_time_step;
3533a2a065fSHong Zhang     break;
3543a2a065fSHong Zhang 
3553a2a065fSHong Zhang   reject_step:
3563a2a065fSHong Zhang     ts->reject++;
3573a2a065fSHong Zhang     accept = PETSC_FALSE;
3583a2a065fSHong Zhang     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
3593a2a065fSHong Zhang       ts->reason = TS_DIVERGED_STEP_REJECTED;
3603a2a065fSHong Zhang       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
3613a2a065fSHong Zhang     }
3623a2a065fSHong Zhang   }
3633a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
3643a2a065fSHong Zhang }
3653a2a065fSHong Zhang 
3663a2a065fSHong Zhang static PetscErrorCode TSSetUp_ARKIMEX_FastSlowSplit(TS ts)
3673a2a065fSHong Zhang {
3683a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
3693a2a065fSHong Zhang   ARKTableau  tab = ark->tableau;
3703a2a065fSHong Zhang   Vec         Xfast, Xslow;
3713a2a065fSHong Zhang 
3723a2a065fSHong Zhang   PetscFunctionBegin;
3733a2a065fSHong Zhang   PetscCall(PetscMalloc1(2 * tab->s, &ark->work));
3743a2a065fSHong Zhang   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y));
3753a2a065fSHong Zhang   PetscCall(TSRHSSplitGetIS(ts, "slow", &ark->is_slow));
3763a2a065fSHong Zhang   PetscCall(TSRHSSplitGetIS(ts, "fast", &ark->is_fast));
3773a2a065fSHong 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");
3783a2a065fSHong Zhang   /* The following vectors need to be resized */
3793a2a065fSHong Zhang   PetscCall(VecDestroy(&ark->Ydot));
3803a2a065fSHong Zhang   PetscCall(VecDestroy(&ark->Ydot0));
3813a2a065fSHong Zhang   PetscCall(VecDestroy(&ark->Z));
3823a2a065fSHong Zhang   PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_fast));
3833a2a065fSHong 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)
3843a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->Y_prev));
3853a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_prev));
3863a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_prev));
3873a2a065fSHong Zhang   }
3883a2a065fSHong Zhang   /* Allocate working vectors */
3893a2a065fSHong Zhang   if (ark->is_fast && ark->is_slow) PetscCall(VecDuplicate(ts->vec_sol, &ark->Y_snes)); // need an additional work vector for SNES
3903a2a065fSHong Zhang   if (ark->is_fast) {
3913a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_fast, &Xfast));
3923a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotRHS_fast));
3933a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotI_fast));
3943a2a065fSHong Zhang     PetscCall(VecDuplicate(Xfast, &ark->Ydot));
3953a2a065fSHong Zhang     PetscCall(VecDuplicate(Xfast, &ark->Ydot0));
3963a2a065fSHong Zhang     PetscCall(VecDuplicate(Xfast, &ark->Z));
3973a2a065fSHong Zhang     if (ark->extrapolate) {
3983a2a065fSHong Zhang       PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->Y_prev));
3993a2a065fSHong Zhang       PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotI_prev));
4003a2a065fSHong Zhang       PetscCall(VecDuplicateVecs(Xfast, tab->s, &ark->YdotRHS_prev));
4013a2a065fSHong Zhang     }
4023a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_fast, &Xfast));
4033a2a065fSHong Zhang   }
4043a2a065fSHong Zhang   if (ark->is_slow) {
4053a2a065fSHong Zhang     PetscCall(VecGetSubVector(ts->vec_sol, ark->is_slow, &Xslow));
4063a2a065fSHong Zhang     PetscCall(VecDuplicateVecs(Xslow, tab->s, &ark->YdotRHS_slow));
4073a2a065fSHong Zhang     PetscCall(VecRestoreSubVector(ts->vec_sol, ark->is_slow, &Xslow));
4083a2a065fSHong Zhang   }
4093a2a065fSHong Zhang   ts->ops->step         = TSStep_ARKIMEX_FastSlowSplit;
4103a2a065fSHong Zhang   ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX_FastSlowSplit;
4113a2a065fSHong Zhang   PetscCall(TSARKIMEXSetSplits(ts));
4123a2a065fSHong 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
4133a2a065fSHong Zhang     SNES snes, snes_fast;
4143a2a065fSHong Zhang     PetscErrorCode (*func)(SNES, Vec, Mat, Mat, void *);
415*4bd3aaa3SHong Zhang     PetscCall(TSRHSSplitGetSNES(ts, &snes));
4163a2a065fSHong Zhang     PetscCall(TSGetSNES(ark->subts_fast, &snes_fast));
4173a2a065fSHong Zhang     PetscCall(SNESGetJacobian(snes_fast, NULL, NULL, &func, NULL));
4183a2a065fSHong Zhang     if (func == SNESTSFormJacobian) PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESTSFormJacobian, ts));
4193a2a065fSHong Zhang   }
4203a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
4213a2a065fSHong Zhang }
4223a2a065fSHong Zhang 
4233a2a065fSHong Zhang static PetscErrorCode TSReset_ARKIMEX_FastSlowSplit(TS ts)
4243a2a065fSHong Zhang {
4253a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
4263a2a065fSHong Zhang   ARKTableau  tab = ark->tableau;
4273a2a065fSHong Zhang 
4283a2a065fSHong Zhang   PetscFunctionBegin;
4293a2a065fSHong Zhang   if (tab) {
4303a2a065fSHong Zhang     PetscCall(PetscFree(ark->work));
4313a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->Y));
4323a2a065fSHong Zhang     if (ark->is_fast && ark->is_slow) PetscCall(VecDestroy(&ark->Y_snes));
4333a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_slow));
4343a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_fast));
4353a2a065fSHong Zhang     PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_fast));
4363a2a065fSHong Zhang     PetscCall(VecDestroy(&ark->Ydot));
4373a2a065fSHong Zhang     PetscCall(VecDestroy(&ark->Ydot0));
4383a2a065fSHong Zhang     PetscCall(VecDestroy(&ark->Z));
4393a2a065fSHong Zhang     if (ark->extrapolate) {
4403a2a065fSHong Zhang       PetscCall(VecDestroyVecs(tab->s, &ark->Y_prev));
4413a2a065fSHong Zhang       PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_prev));
4423a2a065fSHong Zhang       PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_prev));
4433a2a065fSHong Zhang     }
4443a2a065fSHong Zhang   }
4453a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
4463a2a065fSHong Zhang }
4473a2a065fSHong Zhang 
4483a2a065fSHong Zhang PetscErrorCode TSARKIMEXSetFastSlowSplit_ARKIMEX(TS ts, PetscBool fastslowsplit)
4493a2a065fSHong Zhang {
4503a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
4513a2a065fSHong Zhang 
4523a2a065fSHong Zhang   PetscFunctionBegin;
4533a2a065fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4543a2a065fSHong Zhang   ark->fastslowsplit = fastslowsplit;
4553a2a065fSHong Zhang   if (fastslowsplit) {
4563a2a065fSHong Zhang     PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_ARKIMEX_FastSlowSplit_C", TSSetUp_ARKIMEX_FastSlowSplit));
4573a2a065fSHong Zhang     PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_ARKIMEX_FastSlowSplit_C", TSReset_ARKIMEX_FastSlowSplit));
4583a2a065fSHong Zhang   }
4593a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
4603a2a065fSHong Zhang }
4613a2a065fSHong Zhang 
4623a2a065fSHong Zhang PetscErrorCode TSARKIMEXGetFastSlowSplit_ARKIMEX(TS ts, PetscBool *fastslowsplit)
4633a2a065fSHong Zhang {
4643a2a065fSHong Zhang   TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data;
4653a2a065fSHong Zhang 
4663a2a065fSHong Zhang   PetscFunctionBegin;
4673a2a065fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4683a2a065fSHong Zhang   *fastslowsplit = ark->fastslowsplit;
4693a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
4703a2a065fSHong Zhang }
471