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 15abadfa0fSMatthew G. Knepley static PetscErrorCode TSMimexGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 16abadfa0fSMatthew G. Knepley { 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)); 22abadfa0fSMatthew G. Knepley 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)); 26abadfa0fSMatthew G. Knepley else {*Xdot = mimex->Xdot;} 27abadfa0fSMatthew G. Knepley } 28abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 29abadfa0fSMatthew G. Knepley } 30abadfa0fSMatthew G. Knepley 31abadfa0fSMatthew G. Knepley static PetscErrorCode TSMimexRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 32abadfa0fSMatthew G. Knepley { 33abadfa0fSMatthew G. Knepley PetscFunctionBegin; 349566063dSJacob Faibussowitsch if (X0) if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_X0", X0)); 359566063dSJacob Faibussowitsch if (Xdot) if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_Xdot", Xdot)); 36abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 37abadfa0fSMatthew G. Knepley } 38abadfa0fSMatthew G. Knepley 3918c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexGetXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 4018c55e3fSMatthew G. Knepley { 4118c55e3fSMatthew G. Knepley PetscFunctionBegin; 429566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSMimex_Xstar", Xstar)); 439566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSMimex_G", G)); 4418c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 4518c55e3fSMatthew G. Knepley } 4618c55e3fSMatthew G. Knepley 4718c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexRestoreXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 4818c55e3fSMatthew G. Knepley { 4918c55e3fSMatthew G. Knepley PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_Xstar", Xstar)); 519566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSMimex_G", G)); 5218c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 5318c55e3fSMatthew G. Knepley } 5418c55e3fSMatthew G. Knepley 55abadfa0fSMatthew G. Knepley /* 56abadfa0fSMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 57abadfa0fSMatthew G. Knepley G(U) = F[t0+dt, U, (U-U0)*shift] = 0 58abadfa0fSMatthew G. Knepley */ 59abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_Mimex(SNES snes, Vec x, Vec y, TS ts) 60abadfa0fSMatthew G. Knepley { 61abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 62abadfa0fSMatthew G. Knepley DM dm, dmsave; 63abadfa0fSMatthew G. Knepley Vec X0, Xdot; 64abadfa0fSMatthew G. Knepley PetscReal shift = 1./ts->time_step; 65abadfa0fSMatthew G. Knepley 66abadfa0fSMatthew G. Knepley PetscFunctionBegin; 679566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 689566063dSJacob Faibussowitsch PetscCall(TSMimexGetX0AndXdot(ts, dm, &X0, &Xdot)); 699566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); 70abadfa0fSMatthew G. Knepley 71abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 72abadfa0fSMatthew G. Knepley dmsave = ts->dm; 73abadfa0fSMatthew G. Knepley ts->dm = dm; 749566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, mimex->stage_time, x, Xdot, y, PETSC_TRUE)); 7518c55e3fSMatthew G. Knepley if (mimex->version == 1) { 7618c55e3fSMatthew G. Knepley DM dm; 7718c55e3fSMatthew G. Knepley PetscDS prob; 7818c55e3fSMatthew G. Knepley PetscSection s; 797e01ee02SMatthew G. Knepley Vec Xstar = NULL, G = NULL; 8018c55e3fSMatthew G. Knepley const PetscScalar *ax; 8118c55e3fSMatthew G. Knepley PetscScalar *axstar; 8218c55e3fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 8318c55e3fSMatthew G. Knepley 849566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 859566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 869566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 879566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 899566063dSJacob Faibussowitsch PetscCall(TSMimexGetXstarAndG(ts, dm, &Xstar, &G)); 909566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, Xstar)); 919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &ax)); 929566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xstar, &axstar)); 9318c55e3fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 9418c55e3fSMatthew G. Knepley PetscBool implicit; 9518c55e3fSMatthew G. Knepley 969566063dSJacob Faibussowitsch PetscCall(PetscDSGetImplicit(prob, f, &implicit)); 9718c55e3fSMatthew G. Knepley if (!implicit) continue; 9818c55e3fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 9918c55e3fSMatthew G. Knepley PetscScalar *a, *axs; 10018c55e3fSMatthew G. Knepley PetscInt fdof, fcdof, d; 10118c55e3fSMatthew G. Knepley 1029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1039566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &fcdof)); 1049566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, ax, &a)); 1059566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, p, f, axstar, &axs)); 10618c55e3fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) axs[d] = a[d]; 10718c55e3fSMatthew G. Knepley } 10818c55e3fSMatthew G. Knepley } 1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &ax)); 1109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xstar, &axstar)); 1119566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, ts->ptime, Xstar, G)); 1129566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, G)); 1139566063dSJacob Faibussowitsch PetscCall(TSMimexRestoreXstarAndG(ts, dm, &Xstar, &G)); 11418c55e3fSMatthew G. Knepley } 115abadfa0fSMatthew G. Knepley ts->dm = dmsave; 1169566063dSJacob Faibussowitsch PetscCall(TSMimexRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 117abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 118abadfa0fSMatthew G. Knepley } 119abadfa0fSMatthew G. Knepley 120abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_Mimex(SNES snes, Vec x, Mat A, Mat B, TS ts) 121abadfa0fSMatthew G. Knepley { 122abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 123abadfa0fSMatthew G. Knepley DM dm, dmsave; 124abadfa0fSMatthew G. Knepley Vec Xdot; 125abadfa0fSMatthew G. Knepley PetscReal shift = 1./ts->time_step; 126abadfa0fSMatthew G. Knepley 127abadfa0fSMatthew G. Knepley PetscFunctionBegin; 128abadfa0fSMatthew G. Knepley /* th->Xdot has already been computed in SNESTSFormFunction_Mimex (SNES guarantees this) */ 1299566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1309566063dSJacob Faibussowitsch PetscCall(TSMimexGetX0AndXdot(ts, dm, NULL, &Xdot)); 131abadfa0fSMatthew G. Knepley 132abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 133abadfa0fSMatthew G. Knepley dmsave = ts->dm; 134abadfa0fSMatthew G. Knepley ts->dm = dm; 1359566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, mimex->stage_time, x, Xdot, shift, A, B, PETSC_TRUE)); 136abadfa0fSMatthew G. Knepley ts->dm = dmsave; 1379566063dSJacob Faibussowitsch PetscCall(TSMimexRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 138abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 139abadfa0fSMatthew G. Knepley } 140abadfa0fSMatthew G. Knepley 14118c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Split(TS ts) 142abadfa0fSMatthew G. Knepley { 143abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 144abadfa0fSMatthew G. Knepley DM dm; 145abadfa0fSMatthew G. Knepley PetscDS prob; 146abadfa0fSMatthew G. Knepley PetscSection s; 147abadfa0fSMatthew G. Knepley Vec sol = ts->vec_sol, update = mimex->update; 148abadfa0fSMatthew G. Knepley const PetscScalar *aupdate; 149abadfa0fSMatthew G. Knepley PetscScalar *asol, dt = ts->time_step; 150abadfa0fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 151abadfa0fSMatthew G. Knepley 152abadfa0fSMatthew G. Knepley PetscFunctionBegin; 1539566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1549566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 1559566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 1569566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 1579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1589566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime)); 159abadfa0fSMatthew G. Knepley /* Compute implicit update */ 160abadfa0fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 1619566063dSJacob Faibussowitsch PetscCall(VecCopy(sol, update)); 1629566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, update)); 1639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(update, &aupdate)); 1649566063dSJacob Faibussowitsch PetscCall(VecGetArray(sol, &asol)); 165abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 16618c55e3fSMatthew G. Knepley PetscBool implicit; 167abadfa0fSMatthew G. Knepley 1689566063dSJacob Faibussowitsch PetscCall(PetscDSGetImplicit(prob, f, &implicit)); 16918c55e3fSMatthew G. Knepley if (!implicit) continue; 170abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 171abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 172abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 173abadfa0fSMatthew G. Knepley 1749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1759566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &fcdof)); 1769566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au)); 1779566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, p, f, asol, &as)); 178abadfa0fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) as[d] = au[d]; 179abadfa0fSMatthew G. Knepley } 180abadfa0fSMatthew G. Knepley } 1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(update, &aupdate)); 1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(sol, &asol)); 183abadfa0fSMatthew G. Knepley /* Compute explicit update */ 1849566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, ts->ptime, sol, update)); 1859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(update, &aupdate)); 1869566063dSJacob Faibussowitsch PetscCall(VecGetArray(sol, &asol)); 187abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 18818c55e3fSMatthew G. Knepley PetscBool implicit; 189abadfa0fSMatthew G. Knepley 1909566063dSJacob Faibussowitsch PetscCall(PetscDSGetImplicit(prob, f, &implicit)); 19118c55e3fSMatthew G. Knepley if (implicit) continue; 192abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 193abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 194abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 195abadfa0fSMatthew G. Knepley 1969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &fcdof)); 1989566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au)); 1999566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, p, f, asol, &as)); 200abadfa0fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) as[d] += dt*au[d]; 201abadfa0fSMatthew G. Knepley } 202abadfa0fSMatthew G. Knepley } 2039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(update, &aupdate)); 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(sol, &asol)); 2059566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &sol)); 206abadfa0fSMatthew G. Knepley ts->ptime += ts->time_step; 207abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 208abadfa0fSMatthew G. Knepley } 209abadfa0fSMatthew G. Knepley 21018c55e3fSMatthew G. Knepley /* Evalute F at U and G at U0 for explicit fields and U for implicit fields */ 21118c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Implicit(TS ts) 21218c55e3fSMatthew G. Knepley { 21318c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 21418c55e3fSMatthew G. Knepley Vec sol = ts->vec_sol; 21518c55e3fSMatthew G. Knepley Vec update = mimex->update; 21618c55e3fSMatthew G. Knepley 21718c55e3fSMatthew G. Knepley PetscFunctionBegin; 2189566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime)); 21918c55e3fSMatthew G. Knepley /* Compute implicit update */ 22018c55e3fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 22118c55e3fSMatthew G. Knepley ts->ptime += ts->time_step; 2229566063dSJacob Faibussowitsch PetscCall(VecCopy(sol, update)); 2239566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, update)); 2249566063dSJacob Faibussowitsch PetscCall(VecCopy(update, sol)); 2259566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &sol)); 22618c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 22718c55e3fSMatthew G. Knepley } 22818c55e3fSMatthew G. Knepley 22918c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex(TS ts) 23018c55e3fSMatthew G. Knepley { 23118c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 23218c55e3fSMatthew G. Knepley 23318c55e3fSMatthew G. Knepley PetscFunctionBegin; 23418c55e3fSMatthew G. Knepley switch(mimex->version) { 23518c55e3fSMatthew G. Knepley case 0: 2369566063dSJacob Faibussowitsch PetscCall(TSStep_Mimex_Split(ts)); break; 23718c55e3fSMatthew G. Knepley case 1: 2389566063dSJacob Faibussowitsch PetscCall(TSStep_Mimex_Implicit(ts)); break; 23918c55e3fSMatthew G. Knepley default: 24098921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_OUTOFRANGE, "Unknown MIMEX version %d", mimex->version); 24118c55e3fSMatthew G. Knepley } 24218c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 24318c55e3fSMatthew G. Knepley } 24418c55e3fSMatthew G. Knepley 245abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 246abadfa0fSMatthew G. Knepley 247abadfa0fSMatthew G. Knepley static PetscErrorCode TSSetUp_Mimex(TS ts) 248abadfa0fSMatthew G. Knepley { 249abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 250abadfa0fSMatthew G. Knepley 251abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &mimex->update)); 2539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &mimex->Xdot)); 254abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 255abadfa0fSMatthew G. Knepley } 256abadfa0fSMatthew G. Knepley 257abadfa0fSMatthew G. Knepley static PetscErrorCode TSReset_Mimex(TS ts) 258abadfa0fSMatthew G. Knepley { 259abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 260abadfa0fSMatthew G. Knepley 261abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mimex->update)); 2639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mimex->Xdot)); 264abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 265abadfa0fSMatthew G. Knepley } 266abadfa0fSMatthew G. Knepley 267abadfa0fSMatthew G. Knepley static PetscErrorCode TSDestroy_Mimex(TS ts) 268abadfa0fSMatthew G. Knepley { 269abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2709566063dSJacob Faibussowitsch PetscCall(TSReset_Mimex(ts)); 2719566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 272abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 273abadfa0fSMatthew G. Knepley } 274abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 275abadfa0fSMatthew G. Knepley 2764416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Mimex(PetscOptionItems *PetscOptionsObject, TS ts) 277abadfa0fSMatthew G. Knepley { 27818c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 27918c55e3fSMatthew G. Knepley 280abadfa0fSMatthew G. Knepley PetscFunctionBegin; 281*d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "MIMEX ODE solver options"); 28218c55e3fSMatthew G. Knepley { 2839566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_mimex_version", "Algorithm version", "TSMimexSetVersion", mimex->version, &mimex->version, NULL)); 28418c55e3fSMatthew G. Knepley } 285*d0609cedSBarry Smith PetscOptionsHeadEnd(); 286abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 287abadfa0fSMatthew G. Knepley } 288abadfa0fSMatthew G. Knepley 289abadfa0fSMatthew G. Knepley static PetscErrorCode TSView_Mimex(TS ts,PetscViewer viewer) 290abadfa0fSMatthew G. Knepley { 2915bf398f3SMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 2925bf398f3SMatthew G. Knepley PetscBool iascii; 2935bf398f3SMatthew G. Knepley 294abadfa0fSMatthew G. Knepley PetscFunctionBegin; 2959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 2965bf398f3SMatthew G. Knepley if (iascii) { 2979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Version = %D\n", mimex->version)); 2985bf398f3SMatthew G. Knepley } 299abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 300abadfa0fSMatthew G. Knepley } 301abadfa0fSMatthew G. Knepley 302abadfa0fSMatthew G. Knepley static PetscErrorCode TSInterpolate_Mimex(TS ts,PetscReal t,Vec X) 303abadfa0fSMatthew G. Knepley { 304abadfa0fSMatthew G. Knepley PetscReal alpha = (ts->ptime - t)/ts->time_step; 305abadfa0fSMatthew G. Knepley 306abadfa0fSMatthew G. Knepley PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X)); 308abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 309abadfa0fSMatthew G. Knepley } 310abadfa0fSMatthew G. Knepley 311560360afSLisandro Dalcin static PetscErrorCode TSComputeLinearStability_Mimex(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 312abadfa0fSMatthew G. Knepley { 313abadfa0fSMatthew G. Knepley PetscFunctionBegin; 314abadfa0fSMatthew G. Knepley *yr = 1.0 + xr; 315abadfa0fSMatthew G. Knepley *yi = xi; 316abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 317abadfa0fSMatthew G. Knepley } 318abadfa0fSMatthew G. Knepley /* ------------------------------------------------------------ */ 319abadfa0fSMatthew G. Knepley 320abadfa0fSMatthew G. Knepley /*MC 321abadfa0fSMatthew G. Knepley TSMIMEX - ODE solver using the explicit forward Mimex method 322abadfa0fSMatthew G. Knepley 323abadfa0fSMatthew G. Knepley Level: beginner 324abadfa0fSMatthew G. Knepley 325abadfa0fSMatthew G. Knepley .seealso: TSCreate(), TS, TSSetType(), TSBEULER 326abadfa0fSMatthew G. Knepley 327abadfa0fSMatthew G. Knepley M*/ 328abadfa0fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_Mimex(TS ts) 329abadfa0fSMatthew G. Knepley { 330abadfa0fSMatthew G. Knepley TS_Mimex *mimex; 331abadfa0fSMatthew G. Knepley 332abadfa0fSMatthew G. Knepley PetscFunctionBegin; 333abadfa0fSMatthew G. Knepley ts->ops->setup = TSSetUp_Mimex; 334abadfa0fSMatthew G. Knepley ts->ops->step = TSStep_Mimex; 335abadfa0fSMatthew G. Knepley ts->ops->reset = TSReset_Mimex; 336abadfa0fSMatthew G. Knepley ts->ops->destroy = TSDestroy_Mimex; 337abadfa0fSMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_Mimex; 338abadfa0fSMatthew G. Knepley ts->ops->view = TSView_Mimex; 339abadfa0fSMatthew G. Knepley ts->ops->interpolate = TSInterpolate_Mimex; 340abadfa0fSMatthew G. Knepley ts->ops->linearstability = TSComputeLinearStability_Mimex; 341abadfa0fSMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_Mimex; 342abadfa0fSMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_Mimex; 3432ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 344abadfa0fSMatthew G. Knepley 3459566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&mimex)); 346abadfa0fSMatthew G. Knepley ts->data = (void*)mimex; 347fdc0a3ceSMatthew G. Knepley 348fdc0a3ceSMatthew G. Knepley mimex->version = 1; 349abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 350abadfa0fSMatthew G. Knepley } 351