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