1 /* 2 Code for time stepping with the multi-rate Runge-Kutta method 3 4 Notes: 5 1) The general system is written as 6 Udot = F(t,U) for the nonsplit version of multi-rate RK, 7 user should give the indexes for both slow and fast components; 8 2) The general system is written as 9 Usdot = Fs(t,Us,Uf) 10 Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits, 11 user should partioned RHS by themselves and also provide the indexes for both slow and fast components. 12 */ 13 14 #include <petsc/private/tsimpl.h> 15 #include <petscdm.h> 16 #include <../src/ts/impls/explicit/rk/rk.h> 17 #include <../src/ts/impls/explicit/rk/mrk.h> 18 19 static PetscErrorCode TSReset_RK_MultirateNonsplit(TS ts) { 20 TS_RK *rk = (TS_RK *)ts->data; 21 RKTableau tab = rk->tableau; 22 23 PetscFunctionBegin; 24 PetscCall(VecDestroy(&rk->X0)); 25 PetscCall(VecDestroyVecs(tab->s, &rk->YdotRHS_slow)); 26 PetscFunctionReturn(0); 27 } 28 29 static PetscErrorCode TSInterpolate_RK_MultirateNonsplit(TS ts, PetscReal itime, Vec X) { 30 TS_RK *rk = (TS_RK *)ts->data; 31 PetscInt s = rk->tableau->s, p = rk->tableau->p, i, j; 32 PetscReal h = ts->time_step; 33 PetscReal tt, t; 34 PetscScalar *b; 35 const PetscReal *B = rk->tableau->binterp; 36 37 PetscFunctionBegin; 38 PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name); 39 t = (itime - rk->ptime) / h; 40 PetscCall(PetscMalloc1(s, &b)); 41 for (i = 0; i < s; i++) b[i] = 0; 42 for (j = 0, tt = t; j < p; j++, tt *= t) { 43 for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt; 44 } 45 PetscCall(VecCopy(rk->X0, X)); 46 PetscCall(VecMAXPY(X, s, b, rk->YdotRHS_slow)); 47 PetscCall(PetscFree(b)); 48 PetscFunctionReturn(0); 49 } 50 51 static PetscErrorCode TSStepRefine_RK_MultirateNonsplit(TS ts) { 52 TS previousts, subts; 53 TS_RK *rk = (TS_RK *)ts->data; 54 RKTableau tab = rk->tableau; 55 Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 56 Vec vec_fast, subvec_fast; 57 const PetscInt s = tab->s; 58 const PetscReal *A = tab->A, *c = tab->c; 59 PetscScalar *w = rk->work; 60 PetscInt i, j, k; 61 PetscReal t = ts->ptime, h = ts->time_step; 62 63 PetscFunctionBegin; 64 PetscCall(VecDuplicate(ts->vec_sol, &vec_fast)); 65 previousts = rk->subts_current; 66 PetscCall(TSRHSSplitGetSubTS(rk->subts_current, "fast", &subts)); 67 PetscCall(TSRHSSplitGetSubTS(subts, "fast", &subts)); 68 for (k = 0; k < rk->dtratio; k++) { 69 for (i = 0; i < s; i++) { 70 PetscCall(TSInterpolate_RK_MultirateNonsplit(ts, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i])); 71 for (j = 0; j < i; j++) w[j] = h / rk->dtratio * A[i * s + j]; 72 /* 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 */ 73 PetscCall(VecCopy(ts->vec_sol, vec_fast)); 74 PetscCall(VecMAXPY(vec_fast, i, w, YdotRHS)); 75 /* update the fast components in the stage value */ 76 PetscCall(VecGetSubVector(vec_fast, rk->is_fast, &subvec_fast)); 77 PetscCall(VecISCopy(Y[i], rk->is_fast, SCATTER_FORWARD, subvec_fast)); 78 PetscCall(VecRestoreSubVector(vec_fast, rk->is_fast, &subvec_fast)); 79 /* compute the stage RHS */ 80 PetscCall(TSComputeRHSFunction(ts, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i], YdotRHS[i])); 81 } 82 PetscCall(VecCopy(ts->vec_sol, vec_fast)); 83 PetscCall(TSEvaluateStep(ts, tab->order, vec_fast, NULL)); 84 /* update the fast components in the solution */ 85 PetscCall(VecGetSubVector(vec_fast, rk->is_fast, &subvec_fast)); 86 PetscCall(VecISCopy(ts->vec_sol, rk->is_fast, SCATTER_FORWARD, subvec_fast)); 87 PetscCall(VecRestoreSubVector(vec_fast, rk->is_fast, &subvec_fast)); 88 89 if (subts) { 90 Vec *YdotRHS_copy; 91 PetscCall(VecDuplicateVecs(ts->vec_sol, s, &YdotRHS_copy)); 92 rk->subts_current = rk->subts_fast; 93 ts->ptime = t + k * h / rk->dtratio; 94 ts->time_step = h / rk->dtratio; 95 PetscCall(TSRHSSplitGetIS(rk->subts_current, "fast", &rk->is_fast)); 96 for (i = 0; i < s; i++) { 97 PetscCall(VecCopy(rk->YdotRHS_slow[i], YdotRHS_copy[i])); 98 PetscCall(VecCopy(YdotRHS[i], rk->YdotRHS_slow[i])); 99 } 100 101 PetscCall(TSStepRefine_RK_MultirateNonsplit(ts)); 102 103 rk->subts_current = previousts; 104 ts->ptime = t; 105 ts->time_step = h; 106 PetscCall(TSRHSSplitGetIS(previousts, "fast", &rk->is_fast)); 107 for (i = 0; i < s; i++) PetscCall(VecCopy(YdotRHS_copy[i], rk->YdotRHS_slow[i])); 108 PetscCall(VecDestroyVecs(s, &YdotRHS_copy)); 109 } 110 } 111 PetscCall(VecDestroy(&vec_fast)); 112 PetscFunctionReturn(0); 113 } 114 115 static PetscErrorCode TSStep_RK_MultirateNonsplit(TS ts) { 116 TS_RK *rk = (TS_RK *)ts->data; 117 RKTableau tab = rk->tableau; 118 Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS, *YdotRHS_slow = rk->YdotRHS_slow; 119 Vec stage_slow, sol_slow; /* vectors store the slow components */ 120 Vec subvec_slow; /* sub vector to store the slow components */ 121 IS is_slow = rk->is_slow; 122 const PetscInt s = tab->s; 123 const PetscReal *A = tab->A, *c = tab->c; 124 PetscScalar *w = rk->work; 125 PetscInt i, j, dtratio = rk->dtratio; 126 PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step; 127 128 PetscFunctionBegin; 129 rk->status = TS_STEP_INCOMPLETE; 130 PetscCall(VecDuplicate(ts->vec_sol, &stage_slow)); 131 PetscCall(VecDuplicate(ts->vec_sol, &sol_slow)); 132 PetscCall(VecCopy(ts->vec_sol, rk->X0)); 133 for (i = 0; i < s; i++) { 134 rk->stage_time = t + h * c[i]; 135 PetscCall(TSPreStage(ts, rk->stage_time)); 136 PetscCall(VecCopy(ts->vec_sol, Y[i])); 137 for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 138 PetscCall(VecMAXPY(Y[i], i, w, YdotRHS_slow)); 139 PetscCall(TSPostStage(ts, rk->stage_time, i, Y)); 140 /* compute the stage RHS */ 141 PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS_slow[i])); 142 } 143 /* update the slow components in the solution */ 144 rk->YdotRHS = YdotRHS_slow; 145 rk->dtratio = 1; 146 PetscCall(TSEvaluateStep(ts, tab->order, sol_slow, NULL)); 147 rk->dtratio = dtratio; 148 rk->YdotRHS = YdotRHS; 149 /* update the slow components in the solution */ 150 PetscCall(VecGetSubVector(sol_slow, is_slow, &subvec_slow)); 151 PetscCall(VecISCopy(ts->vec_sol, is_slow, SCATTER_FORWARD, subvec_slow)); 152 PetscCall(VecRestoreSubVector(sol_slow, is_slow, &subvec_slow)); 153 154 rk->subts_current = ts; 155 rk->ptime = t; 156 rk->time_step = h; 157 PetscCall(TSStepRefine_RK_MultirateNonsplit(ts)); 158 159 ts->ptime = t + ts->time_step; 160 ts->time_step = next_time_step; 161 rk->status = TS_STEP_COMPLETE; 162 163 /* free memory of work vectors */ 164 PetscCall(VecDestroy(&stage_slow)); 165 PetscCall(VecDestroy(&sol_slow)); 166 PetscFunctionReturn(0); 167 } 168 169 static PetscErrorCode TSSetUp_RK_MultirateNonsplit(TS ts) { 170 TS_RK *rk = (TS_RK *)ts->data; 171 RKTableau tab = rk->tableau; 172 173 PetscFunctionBegin; 174 PetscCall(TSRHSSplitGetIS(ts, "slow", &rk->is_slow)); 175 PetscCall(TSRHSSplitGetIS(ts, "fast", &rk->is_fast)); 176 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"); 177 PetscCall(TSRHSSplitGetSubTS(ts, "slow", &rk->subts_slow)); 178 PetscCall(TSRHSSplitGetSubTS(ts, "fast", &rk->subts_fast)); 179 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"); 180 PetscCall(VecDuplicate(ts->vec_sol, &rk->X0)); 181 PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->YdotRHS_slow)); 182 rk->subts_current = rk->subts_fast; 183 184 ts->ops->step = TSStep_RK_MultirateNonsplit; 185 ts->ops->interpolate = TSInterpolate_RK_MultirateNonsplit; 186 PetscFunctionReturn(0); 187 } 188 189 /* 190 Copy DM from tssrc to tsdest, while keeping the original DMTS and DMSNES in tsdest. 191 */ 192 static PetscErrorCode TSCopyDM(TS tssrc, TS tsdest) { 193 DM newdm, dmsrc, dmdest; 194 195 PetscFunctionBegin; 196 PetscCall(TSGetDM(tssrc, &dmsrc)); 197 PetscCall(DMClone(dmsrc, &newdm)); 198 PetscCall(TSGetDM(tsdest, &dmdest)); 199 PetscCall(DMCopyDMTS(dmdest, newdm)); 200 PetscCall(DMCopyDMSNES(dmdest, newdm)); 201 PetscCall(TSSetDM(tsdest, newdm)); 202 PetscCall(DMDestroy(&newdm)); 203 PetscFunctionReturn(0); 204 } 205 206 static PetscErrorCode TSReset_RK_MultirateSplit(TS ts) { 207 TS_RK *rk = (TS_RK *)ts->data; 208 209 PetscFunctionBegin; 210 if (rk->subts_slow) { 211 PetscCall(PetscFree(rk->subts_slow->data)); 212 rk->subts_slow = NULL; 213 } 214 if (rk->subts_fast) { 215 PetscCall(PetscFree(rk->YdotRHS_fast)); 216 PetscCall(PetscFree(rk->YdotRHS_slow)); 217 PetscCall(VecDestroy(&rk->X0)); 218 PetscCall(TSReset_RK_MultirateSplit(rk->subts_fast)); 219 PetscCall(PetscFree(rk->subts_fast->data)); 220 rk->subts_fast = NULL; 221 } 222 PetscFunctionReturn(0); 223 } 224 225 static PetscErrorCode TSInterpolate_RK_MultirateSplit(TS ts, PetscReal itime, Vec X) { 226 TS_RK *rk = (TS_RK *)ts->data; 227 Vec Xslow; 228 PetscInt s = rk->tableau->s, p = rk->tableau->p, i, j; 229 PetscReal h; 230 PetscReal tt, t; 231 PetscScalar *b; 232 const PetscReal *B = rk->tableau->binterp; 233 234 PetscFunctionBegin; 235 PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name); 236 237 switch (rk->status) { 238 case TS_STEP_INCOMPLETE: 239 case TS_STEP_PENDING: 240 h = ts->time_step; 241 t = (itime - ts->ptime) / h; 242 break; 243 case TS_STEP_COMPLETE: 244 h = ts->ptime - ts->ptime_prev; 245 t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 246 break; 247 default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 248 } 249 PetscCall(PetscMalloc1(s, &b)); 250 for (i = 0; i < s; i++) b[i] = 0; 251 for (j = 0, tt = t; j < p; j++, tt *= t) { 252 for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt; 253 } 254 for (i = 0; i < s; i++) PetscCall(VecGetSubVector(rk->YdotRHS[i], rk->is_slow, &rk->YdotRHS_slow[i])); 255 PetscCall(VecGetSubVector(X, rk->is_slow, &Xslow)); 256 PetscCall(VecISCopy(rk->X0, rk->is_slow, SCATTER_REVERSE, Xslow)); 257 PetscCall(VecMAXPY(Xslow, s, b, rk->YdotRHS_slow)); 258 PetscCall(VecRestoreSubVector(X, rk->is_slow, &Xslow)); 259 for (i = 0; i < s; i++) PetscCall(VecRestoreSubVector(rk->YdotRHS[i], rk->is_slow, &rk->YdotRHS_slow[i])); 260 PetscCall(PetscFree(b)); 261 PetscFunctionReturn(0); 262 } 263 264 /* 265 This is for partitioned RHS multirate RK method 266 The step completion formula is 267 268 x1 = x0 + h b^T YdotRHS 269 270 */ 271 static PetscErrorCode TSEvaluateStep_RK_MultirateSplit(TS ts, PetscInt order, Vec X, PetscBool *done) { 272 TS_RK *rk = (TS_RK *)ts->data; 273 RKTableau tab = rk->tableau; 274 Vec Xslow, Xfast; /* subvectors of X which store slow components and fast components respectively */ 275 PetscScalar *w = rk->work; 276 PetscReal h = ts->time_step; 277 PetscInt s = tab->s, j; 278 279 PetscFunctionBegin; 280 PetscCall(VecCopy(ts->vec_sol, X)); 281 if (rk->slow) { 282 for (j = 0; j < s; j++) w[j] = h * tab->b[j]; 283 PetscCall(VecGetSubVector(ts->vec_sol, rk->is_slow, &Xslow)); 284 PetscCall(VecMAXPY(Xslow, s, w, rk->YdotRHS_slow)); 285 PetscCall(VecRestoreSubVector(ts->vec_sol, rk->is_slow, &Xslow)); 286 } else { 287 for (j = 0; j < s; j++) w[j] = h / rk->dtratio * tab->b[j]; 288 PetscCall(VecGetSubVector(X, rk->is_fast, &Xfast)); 289 PetscCall(VecMAXPY(Xfast, s, w, rk->YdotRHS_fast)); 290 PetscCall(VecRestoreSubVector(X, rk->is_fast, &Xfast)); 291 } 292 PetscFunctionReturn(0); 293 } 294 295 static PetscErrorCode TSStepRefine_RK_MultirateSplit(TS ts) { 296 TS_RK *rk = (TS_RK *)ts->data; 297 TS subts_fast = rk->subts_fast, currentlevelts; 298 TS_RK *subrk_fast = (TS_RK *)subts_fast->data; 299 RKTableau tab = rk->tableau; 300 Vec *Y = rk->Y; 301 Vec *YdotRHS = rk->YdotRHS, *YdotRHS_fast = rk->YdotRHS_fast; 302 Vec Yfast, Xfast; 303 const PetscInt s = tab->s; 304 const PetscReal *A = tab->A, *c = tab->c; 305 PetscScalar *w = rk->work; 306 PetscInt i, j, k; 307 PetscReal t = ts->ptime, h = ts->time_step; 308 309 PetscFunctionBegin; 310 for (k = 0; k < rk->dtratio; k++) { 311 PetscCall(VecGetSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 312 for (i = 0; i < s; i++) PetscCall(VecGetSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 313 /* propagate fast component using small time steps */ 314 for (i = 0; i < s; i++) { 315 /* stage value for slow components */ 316 PetscCall(TSInterpolate_RK_MultirateSplit(rk->ts_root, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i])); 317 currentlevelts = rk->ts_root; 318 while (currentlevelts != ts) { /* all the slow parts need to be interpolated separated */ 319 currentlevelts = ((TS_RK *)currentlevelts->data)->subts_fast; 320 PetscCall(TSInterpolate_RK_MultirateSplit(currentlevelts, t + k * h / rk->dtratio + h / rk->dtratio * c[i], Y[i])); 321 } 322 for (j = 0; j < i; j++) w[j] = h / rk->dtratio * A[i * s + j]; 323 subrk_fast->stage_time = t + h / rk->dtratio * c[i]; 324 PetscCall(TSPreStage(subts_fast, subrk_fast->stage_time)); 325 /* stage value for fast components */ 326 PetscCall(VecGetSubVector(Y[i], rk->is_fast, &Yfast)); 327 PetscCall(VecCopy(Xfast, Yfast)); 328 PetscCall(VecMAXPY(Yfast, i, w, YdotRHS_fast)); 329 PetscCall(VecRestoreSubVector(Y[i], rk->is_fast, &Yfast)); 330 PetscCall(TSPostStage(subts_fast, subrk_fast->stage_time, i, Y)); 331 /* compute the stage RHS for fast components */ 332 PetscCall(TSComputeRHSFunction(subts_fast, t + k * h * rk->dtratio + h / rk->dtratio * c[i], Y[i], YdotRHS_fast[i])); 333 } 334 PetscCall(VecRestoreSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 335 /* update the value of fast components using fast time step */ 336 rk->slow = PETSC_FALSE; 337 PetscCall(TSEvaluateStep_RK_MultirateSplit(ts, tab->order, ts->vec_sol, NULL)); 338 for (i = 0; i < s; i++) PetscCall(VecRestoreSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 339 340 if (subrk_fast->subts_fast) { 341 subts_fast->ptime = t + k * h / rk->dtratio; 342 subts_fast->time_step = h / rk->dtratio; 343 PetscCall(TSStepRefine_RK_MultirateSplit(subts_fast)); 344 } 345 /* update the fast components of the solution */ 346 PetscCall(VecGetSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 347 PetscCall(VecISCopy(rk->X0, rk->is_fast, SCATTER_FORWARD, Xfast)); 348 PetscCall(VecRestoreSubVector(ts->vec_sol, rk->is_fast, &Xfast)); 349 } 350 PetscFunctionReturn(0); 351 } 352 353 static PetscErrorCode TSStep_RK_MultirateSplit(TS ts) { 354 TS_RK *rk = (TS_RK *)ts->data; 355 RKTableau tab = rk->tableau; 356 Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 357 Vec *YdotRHS_fast = rk->YdotRHS_fast, *YdotRHS_slow = rk->YdotRHS_slow; 358 Vec Yslow, Yfast; /* subvectors store the stges of slow components and fast components respectively */ 359 const PetscInt s = tab->s; 360 const PetscReal *A = tab->A, *c = tab->c; 361 PetscScalar *w = rk->work; 362 PetscInt i, j; 363 PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step; 364 365 PetscFunctionBegin; 366 rk->status = TS_STEP_INCOMPLETE; 367 for (i = 0; i < s; i++) { 368 PetscCall(VecGetSubVector(YdotRHS[i], rk->is_slow, &YdotRHS_slow[i])); 369 PetscCall(VecGetSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 370 } 371 PetscCall(VecCopy(ts->vec_sol, rk->X0)); 372 /* propagate both slow and fast components using large time steps */ 373 for (i = 0; i < s; i++) { 374 rk->stage_time = t + h * c[i]; 375 PetscCall(TSPreStage(ts, rk->stage_time)); 376 PetscCall(VecCopy(ts->vec_sol, Y[i])); 377 PetscCall(VecGetSubVector(Y[i], rk->is_fast, &Yfast)); 378 PetscCall(VecGetSubVector(Y[i], rk->is_slow, &Yslow)); 379 for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 380 PetscCall(VecMAXPY(Yfast, i, w, YdotRHS_fast)); 381 PetscCall(VecMAXPY(Yslow, i, w, YdotRHS_slow)); 382 PetscCall(VecRestoreSubVector(Y[i], rk->is_fast, &Yfast)); 383 PetscCall(VecRestoreSubVector(Y[i], rk->is_slow, &Yslow)); 384 PetscCall(TSPostStage(ts, rk->stage_time, i, Y)); 385 PetscCall(TSComputeRHSFunction(rk->subts_slow, t + h * c[i], Y[i], YdotRHS_slow[i])); 386 PetscCall(TSComputeRHSFunction(rk->subts_fast, t + h * c[i], Y[i], YdotRHS_fast[i])); 387 } 388 rk->slow = PETSC_TRUE; 389 /* update the slow components of the solution using slow time step */ 390 PetscCall(TSEvaluateStep_RK_MultirateSplit(ts, tab->order, ts->vec_sol, NULL)); 391 /* YdotRHS will be used for interpolation during refinement */ 392 for (i = 0; i < s; i++) { 393 PetscCall(VecRestoreSubVector(YdotRHS[i], rk->is_slow, &YdotRHS_slow[i])); 394 PetscCall(VecRestoreSubVector(YdotRHS[i], rk->is_fast, &YdotRHS_fast[i])); 395 } 396 397 PetscCall(TSStepRefine_RK_MultirateSplit(ts)); 398 399 ts->ptime = t + ts->time_step; 400 ts->time_step = next_time_step; 401 rk->status = TS_STEP_COMPLETE; 402 PetscFunctionReturn(0); 403 } 404 405 static PetscErrorCode TSSetUp_RK_MultirateSplit(TS ts) { 406 TS_RK *rk = (TS_RK *)ts->data, *nextlevelrk, *currentlevelrk; 407 TS nextlevelts; 408 Vec X0; 409 410 PetscFunctionBegin; 411 PetscCall(TSRHSSplitGetIS(ts, "slow", &rk->is_slow)); 412 PetscCall(TSRHSSplitGetIS(ts, "fast", &rk->is_fast)); 413 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"); 414 PetscCall(TSRHSSplitGetSubTS(ts, "slow", &rk->subts_slow)); 415 PetscCall(TSRHSSplitGetSubTS(ts, "fast", &rk->subts_fast)); 416 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"); 417 418 PetscCall(VecDuplicate(ts->vec_sol, &X0)); 419 /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */ 420 currentlevelrk = rk; 421 while (currentlevelrk->subts_fast) { 422 PetscCall(PetscMalloc1(rk->tableau->s, ¤tlevelrk->YdotRHS_fast)); 423 PetscCall(PetscMalloc1(rk->tableau->s, ¤tlevelrk->YdotRHS_slow)); 424 PetscCall(PetscObjectReference((PetscObject)X0)); 425 currentlevelrk->X0 = X0; 426 currentlevelrk->ts_root = ts; 427 428 /* set up the ts for the slow part */ 429 nextlevelts = currentlevelrk->subts_slow; 430 PetscCall(PetscNewLog(nextlevelts, &nextlevelrk)); 431 nextlevelrk->tableau = rk->tableau; 432 nextlevelrk->work = rk->work; 433 nextlevelrk->Y = rk->Y; 434 nextlevelrk->YdotRHS = rk->YdotRHS; 435 nextlevelts->data = (void *)nextlevelrk; 436 PetscCall(TSCopyDM(ts, nextlevelts)); 437 PetscCall(TSSetSolution(nextlevelts, ts->vec_sol)); 438 439 /* set up the ts for the fast part */ 440 nextlevelts = currentlevelrk->subts_fast; 441 PetscCall(PetscNewLog(nextlevelts, &nextlevelrk)); 442 nextlevelrk->tableau = rk->tableau; 443 nextlevelrk->work = rk->work; 444 nextlevelrk->Y = rk->Y; 445 nextlevelrk->YdotRHS = rk->YdotRHS; 446 nextlevelrk->dtratio = rk->dtratio; 447 PetscCall(TSRHSSplitGetIS(nextlevelts, "slow", &nextlevelrk->is_slow)); 448 PetscCall(TSRHSSplitGetSubTS(nextlevelts, "slow", &nextlevelrk->subts_slow)); 449 PetscCall(TSRHSSplitGetIS(nextlevelts, "fast", &nextlevelrk->is_fast)); 450 PetscCall(TSRHSSplitGetSubTS(nextlevelts, "fast", &nextlevelrk->subts_fast)); 451 nextlevelts->data = (void *)nextlevelrk; 452 PetscCall(TSCopyDM(ts, nextlevelts)); 453 PetscCall(TSSetSolution(nextlevelts, ts->vec_sol)); 454 455 currentlevelrk = nextlevelrk; 456 } 457 PetscCall(VecDestroy(&X0)); 458 459 ts->ops->step = TSStep_RK_MultirateSplit; 460 ts->ops->evaluatestep = TSEvaluateStep_RK_MultirateSplit; 461 ts->ops->interpolate = TSInterpolate_RK_MultirateSplit; 462 PetscFunctionReturn(0); 463 } 464 465 PetscErrorCode TSRKSetMultirate_RK(TS ts, PetscBool use_multirate) { 466 TS_RK *rk = (TS_RK *)ts->data; 467 468 PetscFunctionBegin; 469 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 470 rk->use_multirate = use_multirate; 471 if (use_multirate) { 472 rk->dtratio = 2; 473 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", TSSetUp_RK_MultirateSplit)); 474 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", TSReset_RK_MultirateSplit)); 475 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", TSSetUp_RK_MultirateNonsplit)); 476 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", TSReset_RK_MultirateNonsplit)); 477 } else { 478 rk->dtratio = 0; 479 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", NULL)); 480 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", NULL)); 481 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", NULL)); 482 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", NULL)); 483 } 484 PetscFunctionReturn(0); 485 } 486 487 PetscErrorCode TSRKGetMultirate_RK(TS ts, PetscBool *use_multirate) { 488 TS_RK *rk = (TS_RK *)ts->data; 489 490 PetscFunctionBegin; 491 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 492 *use_multirate = rk->use_multirate; 493 PetscFunctionReturn(0); 494 } 495