1abadfa0fSMatthew G. Knepley /* 2abadfa0fSMatthew G. Knepley Code for Timestepping with my makeshift IMEX. 3abadfa0fSMatthew G. Knepley */ 4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5abadfa0fSMatthew G. Knepley #include <petscds.h> 6ea844a1aSMatthew Knepley #include <petscsection.h> 7abadfa0fSMatthew G. Knepley #include <petscdmplex.h> 8abadfa0fSMatthew G. Knepley 9abadfa0fSMatthew G. Knepley typedef struct { 10010837a9SMatthew G. Knepley Vec Xdot, update; 11abadfa0fSMatthew G. Knepley PetscReal stage_time; 1218c55e3fSMatthew G. Knepley PetscInt version; 13abadfa0fSMatthew G. Knepley } TS_Mimex; 14abadfa0fSMatthew G. Knepley 15d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMimexGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 16d71ae5a4SJacob Faibussowitsch { 17abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 18abadfa0fSMatthew G. Knepley 19abadfa0fSMatthew G. Knepley PetscFunctionBegin; 20abadfa0fSMatthew G. Knepley if (X0) { 219566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSMimex_X0", X0)); 22ad540459SPierre Jolivet else *X0 = ts->vec_sol; 23abadfa0fSMatthew G. Knepley } 24abadfa0fSMatthew G. Knepley if (Xdot) { 259566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSMimex_Xdot", Xdot)); 26ad540459SPierre Jolivet else *Xdot = mimex->Xdot; 27abadfa0fSMatthew G. Knepley } 283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29abadfa0fSMatthew G. Knepley } 30abadfa0fSMatthew G. Knepley 31d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMimexRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 32d71ae5a4SJacob Faibussowitsch { 33abadfa0fSMatthew G. Knepley PetscFunctionBegin; 349371c9d4SSatish Balay if (X0) 359371c9d4SSatish Balay if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_X0", X0)); 369371c9d4SSatish Balay if (Xdot) 379371c9d4SSatish Balay if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_Xdot", Xdot)); 383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39abadfa0fSMatthew G. Knepley } 40abadfa0fSMatthew G. Knepley 41d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMimexGetXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 42d71ae5a4SJacob Faibussowitsch { 4318c55e3fSMatthew G. Knepley PetscFunctionBegin; 449566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSMimex_Xstar", Xstar)); 459566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSMimex_G", G)); 463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4718c55e3fSMatthew G. Knepley } 4818c55e3fSMatthew G. Knepley 49d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMimexRestoreXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 50d71ae5a4SJacob Faibussowitsch { 5118c55e3fSMatthew G. Knepley PetscFunctionBegin; 529566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_Xstar", Xstar)); 539566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_G", G)); 543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5518c55e3fSMatthew G. Knepley } 5618c55e3fSMatthew G. Knepley 57abadfa0fSMatthew G. Knepley /* 58abadfa0fSMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 59abadfa0fSMatthew G. Knepley G(U) = F[t0+dt, U, (U-U0)*shift] = 0 60abadfa0fSMatthew G. Knepley */ 61d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Mimex(SNES snes, Vec x, Vec y, TS ts) 62d71ae5a4SJacob Faibussowitsch { 63abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 64abadfa0fSMatthew G. Knepley DM dm, dmsave; 65abadfa0fSMatthew G. Knepley Vec X0, Xdot; 66abadfa0fSMatthew G. Knepley PetscReal shift = 1. / ts->time_step; 67abadfa0fSMatthew G. Knepley 68abadfa0fSMatthew G. Knepley PetscFunctionBegin; 699566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 709566063dSJacob Faibussowitsch PetscCall(TSMimexGetX0AndXdot(ts, dm, &X0, &Xdot)); 719566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); 72abadfa0fSMatthew G. Knepley 73abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 74abadfa0fSMatthew G. Knepley dmsave = ts->dm; 75abadfa0fSMatthew G. Knepley ts->dm = dm; 769566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, mimex->stage_time, x, Xdot, y, PETSC_TRUE)); 7718c55e3fSMatthew G. Knepley if (mimex->version == 1) { 7818c55e3fSMatthew G. Knepley DM dm; 7918c55e3fSMatthew G. Knepley PetscDS prob; 8018c55e3fSMatthew G. Knepley PetscSection s; 817e01ee02SMatthew G. Knepley Vec Xstar = NULL, G = NULL; 8218c55e3fSMatthew G. Knepley const PetscScalar *ax; 8318c55e3fSMatthew G. Knepley PetscScalar *axstar; 8418c55e3fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 8518c55e3fSMatthew G. Knepley 869566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 879566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 889566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 899566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 919566063dSJacob Faibussowitsch PetscCall(TSMimexGetXstarAndG(ts, dm, &Xstar, &G)); 929566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, Xstar)); 939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &ax)); 949566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xstar, &axstar)); 9518c55e3fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 9618c55e3fSMatthew G. Knepley PetscBool implicit; 9718c55e3fSMatthew G. Knepley 989566063dSJacob Faibussowitsch PetscCall(PetscDSGetImplicit(prob, f, &implicit)); 9918c55e3fSMatthew G. Knepley if (!implicit) continue; 10018c55e3fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 10118c55e3fSMatthew G. Knepley PetscScalar *a, *axs; 10218c55e3fSMatthew G. Knepley PetscInt fdof, fcdof, d; 10318c55e3fSMatthew G. Knepley 1049566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &fcdof)); 1069566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, ax, &a)); 1079566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, p, f, axstar, &axs)); 10818c55e3fSMatthew G. Knepley for (d = 0; d < fdof - fcdof; ++d) axs[d] = a[d]; 10918c55e3fSMatthew G. Knepley } 11018c55e3fSMatthew G. Knepley } 1119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &ax)); 1129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xstar, &axstar)); 1139566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, ts->ptime, Xstar, G)); 1149566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, G)); 1159566063dSJacob Faibussowitsch PetscCall(TSMimexRestoreXstarAndG(ts, dm, &Xstar, &G)); 11618c55e3fSMatthew G. Knepley } 117abadfa0fSMatthew G. Knepley ts->dm = dmsave; 1189566063dSJacob Faibussowitsch PetscCall(TSMimexRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120abadfa0fSMatthew G. Knepley } 121abadfa0fSMatthew G. Knepley 122d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Mimex(SNES snes, Vec x, Mat A, Mat B, TS ts) 123d71ae5a4SJacob Faibussowitsch { 124abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 125abadfa0fSMatthew G. Knepley DM dm, dmsave; 126abadfa0fSMatthew G. Knepley Vec Xdot; 127abadfa0fSMatthew G. Knepley PetscReal shift = 1. / ts->time_step; 128abadfa0fSMatthew G. Knepley 129abadfa0fSMatthew G. Knepley PetscFunctionBegin; 130abadfa0fSMatthew G. Knepley /* th->Xdot has already been computed in SNESTSFormFunction_Mimex (SNES guarantees this) */ 1319566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1329566063dSJacob Faibussowitsch PetscCall(TSMimexGetX0AndXdot(ts, dm, NULL, &Xdot)); 133abadfa0fSMatthew G. Knepley 134abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 135abadfa0fSMatthew G. Knepley dmsave = ts->dm; 136abadfa0fSMatthew G. Knepley ts->dm = dm; 1379566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, mimex->stage_time, x, Xdot, shift, A, B, PETSC_TRUE)); 138abadfa0fSMatthew G. Knepley ts->dm = dmsave; 1399566063dSJacob Faibussowitsch PetscCall(TSMimexRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141abadfa0fSMatthew G. Knepley } 142abadfa0fSMatthew G. Knepley 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Mimex_Split(TS ts) 144d71ae5a4SJacob Faibussowitsch { 145abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 146abadfa0fSMatthew G. Knepley DM dm; 147abadfa0fSMatthew G. Knepley PetscDS prob; 148abadfa0fSMatthew G. Knepley PetscSection s; 149abadfa0fSMatthew G. Knepley Vec sol = ts->vec_sol, update = mimex->update; 150abadfa0fSMatthew G. Knepley const PetscScalar *aupdate; 151abadfa0fSMatthew G. Knepley PetscScalar *asol, dt = ts->time_step; 152abadfa0fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 153abadfa0fSMatthew G. Knepley 154abadfa0fSMatthew G. Knepley PetscFunctionBegin; 1559566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1569566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 1579566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 1589566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 1599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1609566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime)); 161abadfa0fSMatthew G. Knepley /* Compute implicit update */ 162abadfa0fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 1639566063dSJacob Faibussowitsch PetscCall(VecCopy(sol, update)); 1649566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, update)); 1659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(update, &aupdate)); 1669566063dSJacob Faibussowitsch PetscCall(VecGetArray(sol, &asol)); 167abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 16818c55e3fSMatthew G. Knepley PetscBool implicit; 169abadfa0fSMatthew G. Knepley 1709566063dSJacob Faibussowitsch PetscCall(PetscDSGetImplicit(prob, f, &implicit)); 17118c55e3fSMatthew G. Knepley if (!implicit) continue; 172abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 173abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 174abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 175abadfa0fSMatthew G. Knepley 1769566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1779566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &fcdof)); 1789566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au)); 1799566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, p, f, asol, &as)); 180abadfa0fSMatthew G. Knepley for (d = 0; d < fdof - fcdof; ++d) as[d] = au[d]; 181abadfa0fSMatthew G. Knepley } 182abadfa0fSMatthew G. Knepley } 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(update, &aupdate)); 1849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(sol, &asol)); 185abadfa0fSMatthew G. Knepley /* Compute explicit update */ 1869566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, ts->ptime, sol, update)); 1879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(update, &aupdate)); 1889566063dSJacob Faibussowitsch PetscCall(VecGetArray(sol, &asol)); 189abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 19018c55e3fSMatthew G. Knepley PetscBool implicit; 191abadfa0fSMatthew G. Knepley 1929566063dSJacob Faibussowitsch PetscCall(PetscDSGetImplicit(prob, f, &implicit)); 19318c55e3fSMatthew G. Knepley if (implicit) continue; 194abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 195abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 196abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 197abadfa0fSMatthew G. Knepley 1989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1999566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &fcdof)); 2009566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au)); 2019566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, p, f, asol, &as)); 202abadfa0fSMatthew G. Knepley for (d = 0; d < fdof - fcdof; ++d) as[d] += dt * au[d]; 203abadfa0fSMatthew G. Knepley } 204abadfa0fSMatthew G. Knepley } 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(update, &aupdate)); 2069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(sol, &asol)); 2079566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &sol)); 208abadfa0fSMatthew G. Knepley ts->ptime += ts->time_step; 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210abadfa0fSMatthew G. Knepley } 211abadfa0fSMatthew G. Knepley 212d5b43468SJose E. Roman /* Evaluate F at U and G at U0 for explicit fields and U for implicit fields */ 213d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Mimex_Implicit(TS ts) 214d71ae5a4SJacob Faibussowitsch { 21518c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 21618c55e3fSMatthew G. Knepley Vec sol = ts->vec_sol; 21718c55e3fSMatthew G. Knepley Vec update = mimex->update; 21818c55e3fSMatthew G. Knepley 21918c55e3fSMatthew G. Knepley PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime)); 22118c55e3fSMatthew G. Knepley /* Compute implicit update */ 22218c55e3fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 22318c55e3fSMatthew G. Knepley ts->ptime += ts->time_step; 2249566063dSJacob Faibussowitsch PetscCall(VecCopy(sol, update)); 2259566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, update)); 2269566063dSJacob Faibussowitsch PetscCall(VecCopy(update, sol)); 2279566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &sol)); 2283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22918c55e3fSMatthew G. Knepley } 23018c55e3fSMatthew G. Knepley 231d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Mimex(TS ts) 232d71ae5a4SJacob Faibussowitsch { 23318c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 23418c55e3fSMatthew G. Knepley 23518c55e3fSMatthew G. Knepley PetscFunctionBegin; 23618c55e3fSMatthew G. Knepley switch (mimex->version) { 237d71ae5a4SJacob Faibussowitsch case 0: 238d71ae5a4SJacob Faibussowitsch PetscCall(TSStep_Mimex_Split(ts)); 239d71ae5a4SJacob Faibussowitsch break; 240d71ae5a4SJacob Faibussowitsch case 1: 241d71ae5a4SJacob Faibussowitsch PetscCall(TSStep_Mimex_Implicit(ts)); 242d71ae5a4SJacob Faibussowitsch break; 243d71ae5a4SJacob Faibussowitsch default: 244d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Unknown MIMEX version %" PetscInt_FMT, mimex->version); 24518c55e3fSMatthew G. Knepley } 2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24718c55e3fSMatthew G. Knepley } 24818c55e3fSMatthew G. Knepley 249abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 250abadfa0fSMatthew G. Knepley 251d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Mimex(TS ts) 252d71ae5a4SJacob Faibussowitsch { 253abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 254abadfa0fSMatthew G. Knepley 255abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &mimex->update)); 2579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &mimex->Xdot)); 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 259abadfa0fSMatthew G. Knepley } 260abadfa0fSMatthew G. Knepley 261d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Mimex(TS ts) 262d71ae5a4SJacob Faibussowitsch { 263abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 264abadfa0fSMatthew G. Knepley 265abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mimex->update)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mimex->Xdot)); 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269abadfa0fSMatthew G. Knepley } 270abadfa0fSMatthew G. Knepley 271d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Mimex(TS ts) 272d71ae5a4SJacob Faibussowitsch { 273abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2749566063dSJacob Faibussowitsch PetscCall(TSReset_Mimex(ts)); 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277abadfa0fSMatthew G. Knepley } 278abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 279abadfa0fSMatthew G. Knepley 280d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Mimex(TS ts, PetscOptionItems *PetscOptionsObject) 281d71ae5a4SJacob Faibussowitsch { 28218c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 28318c55e3fSMatthew G. Knepley 284abadfa0fSMatthew G. Knepley PetscFunctionBegin; 285d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "MIMEX ODE solver options"); 286d71ae5a4SJacob Faibussowitsch { 287d71ae5a4SJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_mimex_version", "Algorithm version", "TSMimexSetVersion", mimex->version, &mimex->version, NULL)); 288d71ae5a4SJacob Faibussowitsch } 289d0609cedSBarry Smith PetscOptionsHeadEnd(); 2903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 291abadfa0fSMatthew G. Knepley } 292abadfa0fSMatthew G. Knepley 293d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Mimex(TS ts, PetscViewer viewer) 294d71ae5a4SJacob Faibussowitsch { 2955bf398f3SMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *)ts->data; 2965bf398f3SMatthew G. Knepley PetscBool iascii; 2975bf398f3SMatthew G. Knepley 298abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 30048a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Version = %" PetscInt_FMT "\n", mimex->version)); 3013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 302abadfa0fSMatthew G. Knepley } 303abadfa0fSMatthew G. Knepley 304d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Mimex(TS ts, PetscReal t, Vec X) 305d71ae5a4SJacob Faibussowitsch { 306abadfa0fSMatthew G. Knepley PetscReal alpha = (ts->ptime - t) / ts->time_step; 307abadfa0fSMatthew G. Knepley 308abadfa0fSMatthew G. Knepley PetscFunctionBegin; 3099566063dSJacob Faibussowitsch PetscCall(VecAXPBY(ts->vec_sol, 1.0 - alpha, alpha, X)); 3103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 311abadfa0fSMatthew G. Knepley } 312abadfa0fSMatthew G. Knepley 313d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Mimex(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) 314d71ae5a4SJacob Faibussowitsch { 315abadfa0fSMatthew G. Knepley PetscFunctionBegin; 316abadfa0fSMatthew G. Knepley *yr = 1.0 + xr; 317abadfa0fSMatthew G. Knepley *yi = xi; 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319abadfa0fSMatthew G. Knepley } 320abadfa0fSMatthew G. Knepley /* ------------------------------------------------------------ */ 321abadfa0fSMatthew G. Knepley 322abadfa0fSMatthew G. Knepley /*MC 323abadfa0fSMatthew G. Knepley TSMIMEX - ODE solver using the explicit forward Mimex method 324abadfa0fSMatthew G. Knepley 325abadfa0fSMatthew G. Knepley Level: beginner 326abadfa0fSMatthew G. Knepley 327*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER` 328abadfa0fSMatthew G. Knepley M*/ 329d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Mimex(TS ts) 330d71ae5a4SJacob Faibussowitsch { 331abadfa0fSMatthew G. Knepley TS_Mimex *mimex; 332abadfa0fSMatthew G. Knepley 333abadfa0fSMatthew G. Knepley PetscFunctionBegin; 334abadfa0fSMatthew G. Knepley ts->ops->setup = TSSetUp_Mimex; 335abadfa0fSMatthew G. Knepley ts->ops->step = TSStep_Mimex; 336abadfa0fSMatthew G. Knepley ts->ops->reset = TSReset_Mimex; 337abadfa0fSMatthew G. Knepley ts->ops->destroy = TSDestroy_Mimex; 338abadfa0fSMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_Mimex; 339abadfa0fSMatthew G. Knepley ts->ops->view = TSView_Mimex; 340abadfa0fSMatthew G. Knepley ts->ops->interpolate = TSInterpolate_Mimex; 341abadfa0fSMatthew G. Knepley ts->ops->linearstability = TSComputeLinearStability_Mimex; 342abadfa0fSMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_Mimex; 343abadfa0fSMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_Mimex; 3442ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 345abadfa0fSMatthew G. Knepley 3464dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mimex)); 347abadfa0fSMatthew G. Knepley ts->data = (void *)mimex; 348fdc0a3ceSMatthew G. Knepley 349fdc0a3ceSMatthew G. Knepley mimex->version = 1; 3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 351abadfa0fSMatthew G. Knepley } 352