1474dd773SHong Zhang /* 2474dd773SHong Zhang Code for time stepping with the multi-rate Runge-Kutta method 3474dd773SHong Zhang 4474dd773SHong Zhang Notes: 5474dd773SHong Zhang 1) The general system is written as 6474dd773SHong Zhang Udot = F(t,U) for the nonsplit version of multi-rate RK, 7474dd773SHong Zhang user should give the indexes for both slow and fast components; 8474dd773SHong Zhang 2) The general system is written as 9474dd773SHong Zhang Usdot = Fs(t,Us,Uf) 10474dd773SHong Zhang Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits, 11474dd773SHong Zhang user should partioned RHS by themselves and also provide the indexes for both slow and fast components. 12474dd773SHong Zhang */ 13474dd773SHong Zhang 14474dd773SHong Zhang #include <petsc/private/tsimpl.h> 15474dd773SHong Zhang #include <petscdm.h> 16474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h> 1721052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h> 18474dd773SHong Zhang 199371c9d4SSatish Balay static PetscErrorCode TSReset_RK_MultirateNonsplit(TS ts) { 20e5bffa22SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 21e5bffa22SHong Zhang RKTableau tab = rk->tableau; 22e5bffa22SHong Zhang 23e5bffa22SHong Zhang PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->X0)); 259566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &rk->YdotRHS_slow)); 26e5bffa22SHong Zhang PetscFunctionReturn(0); 27e5bffa22SHong Zhang } 28e5bffa22SHong Zhang 299371c9d4SSatish Balay static PetscErrorCode TSInterpolate_RK_MultirateNonsplit(TS ts, PetscReal itime, Vec X) { 30e5bffa22SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 31e5bffa22SHong Zhang PetscInt s = rk->tableau->s, p = rk->tableau->p, i, j; 3263a6f1b4SHong Zhang PetscReal h = ts->time_step; 33e5bffa22SHong Zhang PetscReal tt, t; 34e5bffa22SHong Zhang PetscScalar *b; 35e5bffa22SHong Zhang const PetscReal *B = rk->tableau->binterp; 36e5bffa22SHong Zhang 37e5bffa22SHong Zhang PetscFunctionBegin; 383c633725SBarry Smith PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name); 3963a6f1b4SHong Zhang t = (itime - rk->ptime) / h; 409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &b)); 41e5bffa22SHong Zhang for (i = 0; i < s; i++) b[i] = 0; 42e5bffa22SHong Zhang for (j = 0, tt = t; j < p; j++, tt *= t) { 43ad540459SPierre Jolivet for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt; 44e5bffa22SHong Zhang } 459566063dSJacob Faibussowitsch PetscCall(VecCopy(rk->X0, X)); 469566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, rk->YdotRHS_slow)); 479566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 48e5bffa22SHong Zhang PetscFunctionReturn(0); 49e5bffa22SHong Zhang } 50e5bffa22SHong Zhang 519371c9d4SSatish Balay static PetscErrorCode TSStepRefine_RK_MultirateNonsplit(TS ts) { 5263a6f1b4SHong Zhang TS previousts, subts; 5363a6f1b4SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 5463a6f1b4SHong Zhang RKTableau tab = rk->tableau; 5563a6f1b4SHong Zhang Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 5663a6f1b4SHong Zhang Vec vec_fast, subvec_fast; 5763a6f1b4SHong Zhang const PetscInt s = tab->s; 5863a6f1b4SHong Zhang const PetscReal *A = tab->A, *c = tab->c; 5963a6f1b4SHong Zhang PetscScalar *w = rk->work; 6063a6f1b4SHong Zhang PetscInt i, j, k; 6163a6f1b4SHong Zhang PetscReal t = ts->ptime, h = ts->time_step; 6263a6f1b4SHong Zhang 63362febeeSStefano Zampini PetscFunctionBegin; 649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &vec_fast)); 6563a6f1b4SHong Zhang previousts = rk->subts_current; 669566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(rk->subts_current, "fast", &subts)); 679566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(subts, "fast", &subts)); 6863a6f1b4SHong Zhang for (k = 0; k < rk->dtratio; k++) { 6963a6f1b4SHong Zhang for (i = 0; i < s; i++) { 709566063dSJacob Faibussowitsch PetscCall(TSInterpolate_RK_MultirateNonsplit(ts, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i])); 7163a6f1b4SHong Zhang for (j = 0; j < i; j++) w[j] = h / rk->dtratio * A[i * s + j]; 7263a6f1b4SHong Zhang /* update the fast components in the stage value, the slow components will be ignored, so we do not care the slow part in vec_fast */ 739566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, vec_fast)); 749566063dSJacob Faibussowitsch PetscCall(VecMAXPY(vec_fast, i, w, YdotRHS)); 7563a6f1b4SHong Zhang /* update the fast components in the stage value */ 769566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vec_fast, rk->is_fast, &subvec_fast)); 779566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[i], rk->is_fast, SCATTER_FORWARD, subvec_fast)); 789566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vec_fast, rk->is_fast, &subvec_fast)); 7963a6f1b4SHong Zhang /* compute the stage RHS */ 809566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i], YdotRHS[i])); 8163a6f1b4SHong Zhang } 829566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, vec_fast)); 839566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, vec_fast, NULL)); 8463a6f1b4SHong Zhang /* update the fast components in the solution */ 859566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vec_fast, rk->is_fast, &subvec_fast)); 869566063dSJacob Faibussowitsch PetscCall(VecISCopy(ts->vec_sol, rk->is_fast, SCATTER_FORWARD, subvec_fast)); 879566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vec_fast, rk->is_fast, &subvec_fast)); 8863a6f1b4SHong Zhang 8963a6f1b4SHong Zhang if (subts) { 9063a6f1b4SHong Zhang Vec *YdotRHS_copy; 919566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, s, &YdotRHS_copy)); 9263a6f1b4SHong Zhang rk->subts_current = rk->subts_fast; 9363a6f1b4SHong Zhang ts->ptime = t + k * h / rk->dtratio; 9463a6f1b4SHong Zhang ts->time_step = h / rk->dtratio; 959566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(rk->subts_current, "fast", &rk->is_fast)); 9663a6f1b4SHong Zhang for (i = 0; i < s; i++) { 979566063dSJacob Faibussowitsch PetscCall(VecCopy(rk->YdotRHS_slow[i], YdotRHS_copy[i])); 989566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotRHS[i], rk->YdotRHS_slow[i])); 9963a6f1b4SHong Zhang } 10063a6f1b4SHong Zhang 1019566063dSJacob Faibussowitsch PetscCall(TSStepRefine_RK_MultirateNonsplit(ts)); 10263a6f1b4SHong Zhang 10363a6f1b4SHong Zhang rk->subts_current = previousts; 10463a6f1b4SHong Zhang ts->ptime = t; 10563a6f1b4SHong Zhang ts->time_step = h; 1069566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(previousts, "fast", &rk->is_fast)); 10748a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(VecCopy(YdotRHS_copy[i], rk->YdotRHS_slow[i])); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(s, &YdotRHS_copy)); 10963a6f1b4SHong Zhang } 11063a6f1b4SHong Zhang } 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_fast)); 11263a6f1b4SHong Zhang PetscFunctionReturn(0); 11363a6f1b4SHong Zhang } 11463a6f1b4SHong Zhang 1159371c9d4SSatish Balay static PetscErrorCode TSStep_RK_MultirateNonsplit(TS ts) { 116e5bffa22SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 117e5bffa22SHong Zhang RKTableau tab = rk->tableau; 118e5bffa22SHong Zhang Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS, *YdotRHS_slow = rk->YdotRHS_slow; 119e5bffa22SHong Zhang Vec stage_slow, sol_slow; /* vectors store the slow components */ 120e5bffa22SHong Zhang Vec subvec_slow; /* sub vector to store the slow components */ 121e5bffa22SHong Zhang IS is_slow = rk->is_slow; 122e5bffa22SHong Zhang const PetscInt s = tab->s; 123e5bffa22SHong Zhang const PetscReal *A = tab->A, *c = tab->c; 124e5bffa22SHong Zhang PetscScalar *w = rk->work; 12563a6f1b4SHong Zhang PetscInt i, j, dtratio = rk->dtratio; 126e5bffa22SHong Zhang PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step; 127e5bffa22SHong Zhang 128e5bffa22SHong Zhang PetscFunctionBegin; 129e5bffa22SHong Zhang rk->status = TS_STEP_INCOMPLETE; 1309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &stage_slow)); 1319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &sol_slow)); 1329566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, rk->X0)); 133e5bffa22SHong Zhang for (i = 0; i < s; i++) { 134e5bffa22SHong Zhang rk->stage_time = t + h * c[i]; 1359566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, rk->stage_time)); 1369566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 137e5bffa22SHong Zhang for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 1389566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotRHS_slow)); 1399566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, rk->stage_time, i, Y)); 140e5bffa22SHong Zhang /* compute the stage RHS */ 1419566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS_slow[i])); 142e5bffa22SHong Zhang } 143e5bffa22SHong Zhang /* update the slow components in the solution */ 144e5bffa22SHong Zhang rk->YdotRHS = YdotRHS_slow; 145e5bffa22SHong Zhang rk->dtratio = 1; 1469566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, sol_slow, NULL)); 147e5bffa22SHong Zhang rk->dtratio = dtratio; 148e5bffa22SHong Zhang rk->YdotRHS = YdotRHS; 149e5bffa22SHong Zhang /* update the slow components in the solution */ 1509566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(sol_slow, is_slow, &subvec_slow)); 1519566063dSJacob Faibussowitsch PetscCall(VecISCopy(ts->vec_sol, is_slow, SCATTER_FORWARD, subvec_slow)); 1529566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(sol_slow, is_slow, &subvec_slow)); 153e5bffa22SHong Zhang 15463a6f1b4SHong Zhang rk->subts_current = ts; 15563a6f1b4SHong Zhang rk->ptime = t; 15663a6f1b4SHong Zhang rk->time_step = h; 1579566063dSJacob Faibussowitsch PetscCall(TSStepRefine_RK_MultirateNonsplit(ts)); 15863a6f1b4SHong Zhang 159e5bffa22SHong Zhang ts->ptime = t + ts->time_step; 160e5bffa22SHong Zhang ts->time_step = next_time_step; 161e5bffa22SHong Zhang rk->status = TS_STEP_COMPLETE; 16263a6f1b4SHong Zhang 163e5bffa22SHong Zhang /* free memory of work vectors */ 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&stage_slow)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sol_slow)); 166e5bffa22SHong Zhang PetscFunctionReturn(0); 167e5bffa22SHong Zhang } 168e5bffa22SHong Zhang 1699371c9d4SSatish Balay static PetscErrorCode TSSetUp_RK_MultirateNonsplit(TS ts) { 1700fe4d17eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 1710fe4d17eSHong Zhang RKTableau tab = rk->tableau; 1720fe4d17eSHong Zhang 1730fe4d17eSHong Zhang PetscFunctionBegin; 1749566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "slow", &rk->is_slow)); 1759566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "fast", &rk->is_fast)); 1763c633725SBarry Smith PetscCheck(rk->is_slow && rk->is_fast, 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"); 1779566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "slow", &rk->subts_slow)); 1789566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "fast", &rk->subts_fast)); 1793c633725SBarry Smith PetscCheck(rk->subts_slow && rk->subts_fast, 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"); 1809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &rk->X0)); 1819566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->YdotRHS_slow)); 1820fe4d17eSHong Zhang rk->subts_current = rk->subts_fast; 1830fe4d17eSHong Zhang 18402555c94SHong Zhang ts->ops->step = TSStep_RK_MultirateNonsplit; 18502555c94SHong Zhang ts->ops->interpolate = TSInterpolate_RK_MultirateNonsplit; 1860fe4d17eSHong Zhang PetscFunctionReturn(0); 1870fe4d17eSHong Zhang } 1880fe4d17eSHong Zhang 18963a6f1b4SHong Zhang /* 19063a6f1b4SHong Zhang Copy DM from tssrc to tsdest, while keeping the original DMTS and DMSNES in tsdest. 19163a6f1b4SHong Zhang */ 1929371c9d4SSatish Balay static PetscErrorCode TSCopyDM(TS tssrc, TS tsdest) { 19363a6f1b4SHong Zhang DM newdm, dmsrc, dmdest; 19463a6f1b4SHong Zhang 19563a6f1b4SHong Zhang PetscFunctionBegin; 1969566063dSJacob Faibussowitsch PetscCall(TSGetDM(tssrc, &dmsrc)); 1979566063dSJacob Faibussowitsch PetscCall(DMClone(dmsrc, &newdm)); 1989566063dSJacob Faibussowitsch PetscCall(TSGetDM(tsdest, &dmdest)); 1999566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(dmdest, newdm)); 2009566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(dmdest, newdm)); 2019566063dSJacob Faibussowitsch PetscCall(TSSetDM(tsdest, newdm)); 2029566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 20363a6f1b4SHong Zhang PetscFunctionReturn(0); 20463a6f1b4SHong Zhang } 20563a6f1b4SHong Zhang 2069371c9d4SSatish Balay static PetscErrorCode TSReset_RK_MultirateSplit(TS ts) { 207474dd773SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 208474dd773SHong Zhang 209474dd773SHong Zhang PetscFunctionBegin; 210bb42530cSHong Zhang if (rk->subts_slow) { 2119566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->subts_slow->data)); 212bb42530cSHong Zhang rk->subts_slow = NULL; 213bb42530cSHong Zhang } 214bb42530cSHong Zhang if (rk->subts_fast) { 2159566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->YdotRHS_fast)); 2169566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->YdotRHS_slow)); 2179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->X0)); 2189566063dSJacob Faibussowitsch PetscCall(TSReset_RK_MultirateSplit(rk->subts_fast)); 2199566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->subts_fast->data)); 220bb42530cSHong Zhang rk->subts_fast = NULL; 221bb42530cSHong Zhang } 222474dd773SHong Zhang PetscFunctionReturn(0); 223474dd773SHong Zhang } 224474dd773SHong Zhang 2259371c9d4SSatish Balay static PetscErrorCode TSInterpolate_RK_MultirateSplit(TS ts, PetscReal itime, Vec X) { 226474dd773SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 22763a6f1b4SHong Zhang Vec Xslow; 228474dd773SHong Zhang PetscInt s = rk->tableau->s, p = rk->tableau->p, i, j; 229474dd773SHong Zhang PetscReal h; 230474dd773SHong Zhang PetscReal tt, t; 231474dd773SHong Zhang PetscScalar *b; 232474dd773SHong Zhang const PetscReal *B = rk->tableau->binterp; 233474dd773SHong Zhang 234474dd773SHong Zhang PetscFunctionBegin; 2353c633725SBarry Smith PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name); 236474dd773SHong Zhang 237474dd773SHong Zhang switch (rk->status) { 238474dd773SHong Zhang case TS_STEP_INCOMPLETE: 239474dd773SHong Zhang case TS_STEP_PENDING: 240474dd773SHong Zhang h = ts->time_step; 241474dd773SHong Zhang t = (itime - ts->ptime) / h; 242474dd773SHong Zhang break; 243474dd773SHong Zhang case TS_STEP_COMPLETE: 244474dd773SHong Zhang h = ts->ptime - ts->ptime_prev; 245474dd773SHong Zhang t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 246474dd773SHong Zhang break; 247474dd773SHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 248474dd773SHong Zhang } 2499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &b)); 250474dd773SHong Zhang for (i = 0; i < s; i++) b[i] = 0; 251474dd773SHong Zhang for (j = 0, tt = t; j < p; j++, tt *= t) { 252ad540459SPierre Jolivet for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt; 253474dd773SHong Zhang } 25448a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(VecGetSubVector(rk->YdotRHS[i], rk->is_slow, &rk->YdotRHS_slow[i])); 2559566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, rk->is_slow, &Xslow)); 2569566063dSJacob Faibussowitsch PetscCall(VecISCopy(rk->X0, rk->is_slow, SCATTER_REVERSE, Xslow)); 2579566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslow, s, b, rk->YdotRHS_slow)); 2589566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, rk->is_slow, &Xslow)); 25948a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(VecRestoreSubVector(rk->YdotRHS[i], rk->is_slow, &rk->YdotRHS_slow[i])); 2609566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 261474dd773SHong Zhang PetscFunctionReturn(0); 262474dd773SHong Zhang } 263474dd773SHong Zhang 264474dd773SHong Zhang /* 265474dd773SHong Zhang This is for partitioned RHS multirate RK method 266474dd773SHong Zhang The step completion formula is 267474dd773SHong Zhang 268474dd773SHong Zhang x1 = x0 + h b^T YdotRHS 269474dd773SHong Zhang 270474dd773SHong Zhang */ 2719371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_RK_MultirateSplit(TS ts, PetscInt order, Vec X, PetscBool *done) { 272474dd773SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 273474dd773SHong Zhang RKTableau tab = rk->tableau; 274474dd773SHong Zhang Vec Xslow, Xfast; /* subvectors of X which store slow components and fast components respectively */ 275474dd773SHong Zhang PetscScalar *w = rk->work; 276474dd773SHong Zhang PetscReal h = ts->time_step; 277474dd773SHong Zhang PetscInt s = tab->s, j; 278474dd773SHong Zhang 279474dd773SHong Zhang PetscFunctionBegin; 2809566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 281474dd773SHong Zhang if (rk->slow) { 282474dd773SHong Zhang for (j = 0; j < s; j++) w[j] = h * tab->b[j]; 2839566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, rk->is_slow, &Xslow)); 2849566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslow, s, w, rk->YdotRHS_slow)); 2859566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, rk->is_slow, &Xslow)); 286474dd773SHong Zhang } else { 287474dd773SHong Zhang for (j = 0; j < s; j++) w[j] = h / rk->dtratio * tab->b[j]; 2889566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, rk->is_fast, &Xfast)); 2899566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xfast, s, w, rk->YdotRHS_fast)); 2909566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, rk->is_fast, &Xfast)); 291474dd773SHong Zhang } 292474dd773SHong Zhang PetscFunctionReturn(0); 293474dd773SHong Zhang } 294474dd773SHong Zhang 2959371c9d4SSatish Balay static PetscErrorCode TSStepRefine_RK_MultirateSplit(TS ts) { 29663a6f1b4SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 29763a6f1b4SHong Zhang TS subts_fast = rk->subts_fast, currentlevelts; 29863a6f1b4SHong Zhang TS_RK *subrk_fast = (TS_RK *)subts_fast->data; 29963a6f1b4SHong Zhang RKTableau tab = rk->tableau; 30063a6f1b4SHong Zhang Vec *Y = rk->Y; 30163a6f1b4SHong Zhang Vec *YdotRHS = rk->YdotRHS, *YdotRHS_fast = rk->YdotRHS_fast; 30263a6f1b4SHong Zhang Vec Yfast, Xfast; 30363a6f1b4SHong Zhang const PetscInt s = tab->s; 30463a6f1b4SHong Zhang const PetscReal *A = tab->A, *c = tab->c; 30563a6f1b4SHong Zhang PetscScalar *w = rk->work; 30663a6f1b4SHong Zhang PetscInt i, j, k; 30763a6f1b4SHong Zhang PetscReal t = ts->ptime, h = ts->time_step; 30863a6f1b4SHong Zhang 309362febeeSStefano Zampini PetscFunctionBegin; 31063a6f1b4SHong Zhang for (k = 0; k < rk->dtratio; k++) { 3119566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 31248a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(VecGetSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 313a5b23f4aSJose E. Roman /* propagate fast component using small time steps */ 31463a6f1b4SHong Zhang for (i = 0; i < s; i++) { 31563a6f1b4SHong Zhang /* stage value for slow components */ 3169566063dSJacob Faibussowitsch PetscCall(TSInterpolate_RK_MultirateSplit(rk->ts_root, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i])); 31763a6f1b4SHong Zhang currentlevelts = rk->ts_root; 31863a6f1b4SHong Zhang while (currentlevelts != ts) { /* all the slow parts need to be interpolated separated */ 31963a6f1b4SHong Zhang currentlevelts = ((TS_RK *)currentlevelts->data)->subts_fast; 3209566063dSJacob Faibussowitsch PetscCall(TSInterpolate_RK_MultirateSplit(currentlevelts, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i])); 32163a6f1b4SHong Zhang } 32263a6f1b4SHong Zhang for (j = 0; j < i; j++) w[j] = h / rk->dtratio * A[i * s + j]; 32363a6f1b4SHong Zhang subrk_fast->stage_time = t + h / rk->dtratio * c[i]; 3249566063dSJacob Faibussowitsch PetscCall(TSPreStage(subts_fast, subrk_fast->stage_time)); 32563a6f1b4SHong Zhang /* stage value for fast components */ 3269566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], rk->is_fast, &Yfast)); 3279566063dSJacob Faibussowitsch PetscCall(VecCopy(Xfast, Yfast)); 3289566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast, i, w, YdotRHS_fast)); 3299566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], rk->is_fast, &Yfast)); 3309566063dSJacob Faibussowitsch PetscCall(TSPostStage(subts_fast, subrk_fast->stage_time, i, Y)); 33163a6f1b4SHong Zhang /* compute the stage RHS for fast components */ 3329566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(subts_fast, t + k * h * rk->dtratio + h / rk->dtratio * c[i], Y[i], YdotRHS_fast[i])); 33363a6f1b4SHong Zhang } 3349566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 33563a6f1b4SHong Zhang /* update the value of fast components using fast time step */ 33663a6f1b4SHong Zhang rk->slow = PETSC_FALSE; 3379566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_RK_MultirateSplit(ts, tab->order, ts->vec_sol, NULL)); 33848a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(VecRestoreSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 33963a6f1b4SHong Zhang 34063a6f1b4SHong Zhang if (subrk_fast->subts_fast) { 34163a6f1b4SHong Zhang subts_fast->ptime = t + k * h / rk->dtratio; 34263a6f1b4SHong Zhang subts_fast->time_step = h / rk->dtratio; 3439566063dSJacob Faibussowitsch PetscCall(TSStepRefine_RK_MultirateSplit(subts_fast)); 34463a6f1b4SHong Zhang } 34563a6f1b4SHong Zhang /* update the fast components of the solution */ 3469566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 3479566063dSJacob Faibussowitsch PetscCall(VecISCopy(rk->X0, rk->is_fast, SCATTER_FORWARD, Xfast)); 3489566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 34963a6f1b4SHong Zhang } 35063a6f1b4SHong Zhang PetscFunctionReturn(0); 35163a6f1b4SHong Zhang } 35263a6f1b4SHong Zhang 3539371c9d4SSatish Balay static PetscErrorCode TSStep_RK_MultirateSplit(TS ts) { 354474dd773SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 355474dd773SHong Zhang RKTableau tab = rk->tableau; 356474dd773SHong Zhang Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 35763a6f1b4SHong Zhang Vec *YdotRHS_fast = rk->YdotRHS_fast, *YdotRHS_slow = rk->YdotRHS_slow; 358474dd773SHong Zhang Vec Yslow, Yfast; /* subvectors store the stges of slow components and fast components respectively */ 359474dd773SHong Zhang const PetscInt s = tab->s; 360474dd773SHong Zhang const PetscReal *A = tab->A, *c = tab->c; 361474dd773SHong Zhang PetscScalar *w = rk->work; 36263a6f1b4SHong Zhang PetscInt i, j; 363474dd773SHong Zhang PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step; 364474dd773SHong Zhang 365474dd773SHong Zhang PetscFunctionBegin; 366474dd773SHong Zhang rk->status = TS_STEP_INCOMPLETE; 367474dd773SHong Zhang for (i = 0; i < s; i++) { 3689566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(YdotRHS[i], rk->is_slow, &YdotRHS_slow[i])); 3699566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 370474dd773SHong Zhang } 3719566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, rk->X0)); 372a5b23f4aSJose E. Roman /* propagate both slow and fast components using large time steps */ 373474dd773SHong Zhang for (i = 0; i < s; i++) { 374474dd773SHong Zhang rk->stage_time = t + h * c[i]; 3759566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, rk->stage_time)); 3769566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 3779566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], rk->is_fast, &Yfast)); 3789566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], rk->is_slow, &Yslow)); 379474dd773SHong Zhang for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 3809566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast, i, w, YdotRHS_fast)); 3819566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslow, i, w, YdotRHS_slow)); 3829566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], rk->is_fast, &Yfast)); 3839566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], rk->is_slow, &Yslow)); 3849566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, rk->stage_time, i, Y)); 3859566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(rk->subts_slow, t + h * c[i], Y[i], YdotRHS_slow[i])); 3869566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(rk->subts_fast, t + h * c[i], Y[i], YdotRHS_fast[i])); 387474dd773SHong Zhang } 388474dd773SHong Zhang rk->slow = PETSC_TRUE; 38963a6f1b4SHong Zhang /* update the slow components of the solution using slow time step */ 3909566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_RK_MultirateSplit(ts, tab->order, ts->vec_sol, NULL)); 39163a6f1b4SHong Zhang /* YdotRHS will be used for interpolation during refinement */ 392474dd773SHong Zhang for (i = 0; i < s; i++) { 3939566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(YdotRHS[i], rk->is_slow, &YdotRHS_slow[i])); 3949566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 395474dd773SHong Zhang } 39663a6f1b4SHong Zhang 3979566063dSJacob Faibussowitsch PetscCall(TSStepRefine_RK_MultirateSplit(ts)); 39863a6f1b4SHong Zhang 399bb42530cSHong Zhang ts->ptime = t + ts->time_step; 400474dd773SHong Zhang ts->time_step = next_time_step; 401474dd773SHong Zhang rk->status = TS_STEP_COMPLETE; 402474dd773SHong Zhang PetscFunctionReturn(0); 403474dd773SHong Zhang } 404474dd773SHong Zhang 4059371c9d4SSatish Balay static PetscErrorCode TSSetUp_RK_MultirateSplit(TS ts) { 4060fe4d17eSHong Zhang TS_RK *rk = (TS_RK *)ts->data, *nextlevelrk, *currentlevelrk; 4070fe4d17eSHong Zhang TS nextlevelts; 4080fe4d17eSHong Zhang Vec X0; 4090fe4d17eSHong Zhang 4100fe4d17eSHong Zhang PetscFunctionBegin; 4119566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "slow", &rk->is_slow)); 4129566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "fast", &rk->is_fast)); 4133c633725SBarry Smith PetscCheck(rk->is_slow && rk->is_fast, 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"); 4149566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "slow", &rk->subts_slow)); 4159566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "fast", &rk->subts_fast)); 4163c633725SBarry Smith PetscCheck(rk->subts_slow && rk->subts_fast, 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"); 4170fe4d17eSHong Zhang 4189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &X0)); 4190fe4d17eSHong Zhang /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */ 4200fe4d17eSHong Zhang currentlevelrk = rk; 4210fe4d17eSHong Zhang while (currentlevelrk->subts_fast) { 4229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rk->tableau->s, ¤tlevelrk->YdotRHS_fast)); 4239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rk->tableau->s, ¤tlevelrk->YdotRHS_slow)); 4249566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)X0)); 4250fe4d17eSHong Zhang currentlevelrk->X0 = X0; 4260fe4d17eSHong Zhang currentlevelrk->ts_root = ts; 4270fe4d17eSHong Zhang 4280fe4d17eSHong Zhang /* set up the ts for the slow part */ 4290fe4d17eSHong Zhang nextlevelts = currentlevelrk->subts_slow; 430*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&nextlevelrk)); 4310fe4d17eSHong Zhang nextlevelrk->tableau = rk->tableau; 4320fe4d17eSHong Zhang nextlevelrk->work = rk->work; 4330fe4d17eSHong Zhang nextlevelrk->Y = rk->Y; 4340fe4d17eSHong Zhang nextlevelrk->YdotRHS = rk->YdotRHS; 4350fe4d17eSHong Zhang nextlevelts->data = (void *)nextlevelrk; 4369566063dSJacob Faibussowitsch PetscCall(TSCopyDM(ts, nextlevelts)); 4379566063dSJacob Faibussowitsch PetscCall(TSSetSolution(nextlevelts, ts->vec_sol)); 4380fe4d17eSHong Zhang 4390fe4d17eSHong Zhang /* set up the ts for the fast part */ 4400fe4d17eSHong Zhang nextlevelts = currentlevelrk->subts_fast; 441*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&nextlevelrk)); 4420fe4d17eSHong Zhang nextlevelrk->tableau = rk->tableau; 4430fe4d17eSHong Zhang nextlevelrk->work = rk->work; 4440fe4d17eSHong Zhang nextlevelrk->Y = rk->Y; 4450fe4d17eSHong Zhang nextlevelrk->YdotRHS = rk->YdotRHS; 4460fe4d17eSHong Zhang nextlevelrk->dtratio = rk->dtratio; 4479566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(nextlevelts, "slow", &nextlevelrk->is_slow)); 4489566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(nextlevelts, "slow", &nextlevelrk->subts_slow)); 4499566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(nextlevelts, "fast", &nextlevelrk->is_fast)); 4509566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(nextlevelts, "fast", &nextlevelrk->subts_fast)); 4510fe4d17eSHong Zhang nextlevelts->data = (void *)nextlevelrk; 4529566063dSJacob Faibussowitsch PetscCall(TSCopyDM(ts, nextlevelts)); 4539566063dSJacob Faibussowitsch PetscCall(TSSetSolution(nextlevelts, ts->vec_sol)); 4540fe4d17eSHong Zhang 4550fe4d17eSHong Zhang currentlevelrk = nextlevelrk; 4560fe4d17eSHong Zhang } 4579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 4580fe4d17eSHong Zhang 45902555c94SHong Zhang ts->ops->step = TSStep_RK_MultirateSplit; 46002555c94SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_RK_MultirateSplit; 46102555c94SHong Zhang ts->ops->interpolate = TSInterpolate_RK_MultirateSplit; 4620fe4d17eSHong Zhang PetscFunctionReturn(0); 4630fe4d17eSHong Zhang } 4640fe4d17eSHong Zhang 4659371c9d4SSatish Balay PetscErrorCode TSRKSetMultirate_RK(TS ts, PetscBool use_multirate) { 466474dd773SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 467474dd773SHong Zhang 468474dd773SHong Zhang PetscFunctionBegin; 469474dd773SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 4700fe4d17eSHong Zhang rk->use_multirate = use_multirate; 4710fe4d17eSHong Zhang if (use_multirate) { 472474dd773SHong Zhang rk->dtratio = 2; 4739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", TSSetUp_RK_MultirateSplit)); 4749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", TSReset_RK_MultirateSplit)); 4759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", TSSetUp_RK_MultirateNonsplit)); 4769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", TSReset_RK_MultirateNonsplit)); 4770fe4d17eSHong Zhang } else { 4780fe4d17eSHong Zhang rk->dtratio = 0; 4799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", NULL)); 4809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", NULL)); 4819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", NULL)); 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", NULL)); 483474dd773SHong Zhang } 484474dd773SHong Zhang PetscFunctionReturn(0); 485474dd773SHong Zhang } 486474dd773SHong Zhang 4879371c9d4SSatish Balay PetscErrorCode TSRKGetMultirate_RK(TS ts, PetscBool *use_multirate) { 4880fe4d17eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 4890fe4d17eSHong Zhang 4900fe4d17eSHong Zhang PetscFunctionBegin; 4910fe4d17eSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 4920fe4d17eSHong Zhang *use_multirate = rk->use_multirate; 4930fe4d17eSHong Zhang PetscFunctionReturn(0); 4940fe4d17eSHong Zhang } 495