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> 6abadfa0fSMatthew G. Knepley #include <petscdmplex.h> 7abadfa0fSMatthew G. Knepley 8abadfa0fSMatthew G. Knepley typedef struct { 9010837a9SMatthew G. Knepley Vec Xdot, update; 10abadfa0fSMatthew G. Knepley PetscReal stage_time; 1118c55e3fSMatthew G. Knepley PetscInt version; 12abadfa0fSMatthew G. Knepley } TS_Mimex; 13abadfa0fSMatthew G. Knepley 14abadfa0fSMatthew G. Knepley static PetscErrorCode TSMimexGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 15abadfa0fSMatthew G. Knepley { 16abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 17abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 18abadfa0fSMatthew G. Knepley 19abadfa0fSMatthew G. Knepley PetscFunctionBegin; 20abadfa0fSMatthew G. Knepley if (X0) { 21abadfa0fSMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSMimex_X0", X0);CHKERRQ(ierr);} 22abadfa0fSMatthew G. Knepley else {*X0 = ts->vec_sol;} 23abadfa0fSMatthew G. Knepley } 24abadfa0fSMatthew G. Knepley if (Xdot) { 25abadfa0fSMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSMimex_Xdot", Xdot);CHKERRQ(ierr);} 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 PetscErrorCode ierr; 34abadfa0fSMatthew G. Knepley 35abadfa0fSMatthew G. Knepley PetscFunctionBegin; 36abadfa0fSMatthew G. Knepley if (X0) if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_X0", X0);CHKERRQ(ierr);} 37abadfa0fSMatthew G. Knepley if (Xdot) if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_Xdot", Xdot);CHKERRQ(ierr);} 38abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 39abadfa0fSMatthew G. Knepley } 40abadfa0fSMatthew G. Knepley 4118c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexGetXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 4218c55e3fSMatthew G. Knepley { 4318c55e3fSMatthew G. Knepley PetscErrorCode ierr; 4418c55e3fSMatthew G. Knepley 4518c55e3fSMatthew G. Knepley PetscFunctionBegin; 46010837a9SMatthew G. Knepley ierr = DMGetNamedGlobalVector(dm, "TSMimex_Xstar", Xstar);CHKERRQ(ierr); 47010837a9SMatthew G. Knepley ierr = DMGetNamedGlobalVector(dm, "TSMimex_G", G);CHKERRQ(ierr); 4818c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 4918c55e3fSMatthew G. Knepley } 5018c55e3fSMatthew G. Knepley 5118c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexRestoreXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 5218c55e3fSMatthew G. Knepley { 53010837a9SMatthew G. Knepley PetscErrorCode ierr; 54010837a9SMatthew G. Knepley 5518c55e3fSMatthew G. Knepley PetscFunctionBegin; 56010837a9SMatthew G. Knepley ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_Xstar", Xstar);CHKERRQ(ierr); 57010837a9SMatthew G. Knepley ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_G", G);CHKERRQ(ierr); 5818c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 5918c55e3fSMatthew G. Knepley } 6018c55e3fSMatthew G. Knepley 61abadfa0fSMatthew G. Knepley /* 62abadfa0fSMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 63abadfa0fSMatthew G. Knepley G(U) = F[t0+dt, U, (U-U0)*shift] = 0 64abadfa0fSMatthew G. Knepley */ 65abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_Mimex(SNES snes, Vec x, Vec y, TS ts) 66abadfa0fSMatthew G. Knepley { 67abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 68abadfa0fSMatthew G. Knepley DM dm, dmsave; 69abadfa0fSMatthew G. Knepley Vec X0, Xdot; 70abadfa0fSMatthew G. Knepley PetscReal shift = 1./ts->time_step; 71abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 72abadfa0fSMatthew G. Knepley 73abadfa0fSMatthew G. Knepley PetscFunctionBegin; 74abadfa0fSMatthew G. Knepley ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 75abadfa0fSMatthew G. Knepley ierr = TSMimexGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 76abadfa0fSMatthew G. Knepley ierr = VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x);CHKERRQ(ierr); 77abadfa0fSMatthew G. Knepley 78abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 79abadfa0fSMatthew G. Knepley dmsave = ts->dm; 80abadfa0fSMatthew G. Knepley ts->dm = dm; 81abadfa0fSMatthew G. Knepley ierr = TSComputeIFunction(ts, mimex->stage_time, x, Xdot, y, PETSC_TRUE);CHKERRQ(ierr); 8218c55e3fSMatthew G. Knepley if (mimex->version == 1) { 8318c55e3fSMatthew G. Knepley DM dm; 8418c55e3fSMatthew G. Knepley PetscDS prob; 8518c55e3fSMatthew G. Knepley PetscSection s; 867e01ee02SMatthew G. Knepley Vec Xstar = NULL, G = NULL; 8718c55e3fSMatthew G. Knepley const PetscScalar *ax; 8818c55e3fSMatthew G. Knepley PetscScalar *axstar; 8918c55e3fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 9018c55e3fSMatthew G. Knepley 9118c55e3fSMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 9218c55e3fSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 93*92fd8e1eSJed Brown ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 9418c55e3fSMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 9518c55e3fSMatthew G. Knepley ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 9618c55e3fSMatthew G. Knepley ierr = TSMimexGetXstarAndG(ts, dm, &Xstar, &G);CHKERRQ(ierr); 9718c55e3fSMatthew G. Knepley ierr = VecCopy(X0, Xstar);CHKERRQ(ierr); 9818c55e3fSMatthew G. Knepley ierr = VecGetArrayRead(x, &ax);CHKERRQ(ierr); 9918c55e3fSMatthew G. Knepley ierr = VecGetArray(Xstar, &axstar);CHKERRQ(ierr); 10018c55e3fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 10118c55e3fSMatthew G. Knepley PetscBool implicit; 10218c55e3fSMatthew G. Knepley 10318c55e3fSMatthew G. Knepley ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr); 10418c55e3fSMatthew G. Knepley if (!implicit) continue; 10518c55e3fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 10618c55e3fSMatthew G. Knepley PetscScalar *a, *axs; 10718c55e3fSMatthew G. Knepley PetscInt fdof, fcdof, d; 10818c55e3fSMatthew G. Knepley 10918c55e3fSMatthew G. Knepley ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 11018c55e3fSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr); 11118c55e3fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRead(dm, p, f, ax, &a);CHKERRQ(ierr); 11218c55e3fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRef(dm, p, f, axstar, &axs);CHKERRQ(ierr); 11318c55e3fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) axs[d] = a[d]; 11418c55e3fSMatthew G. Knepley } 11518c55e3fSMatthew G. Knepley } 11618c55e3fSMatthew G. Knepley ierr = VecRestoreArrayRead(x, &ax);CHKERRQ(ierr); 11718c55e3fSMatthew G. Knepley ierr = VecRestoreArray(Xstar, &axstar);CHKERRQ(ierr); 11818c55e3fSMatthew G. Knepley ierr = TSComputeRHSFunction(ts, ts->ptime, Xstar, G);CHKERRQ(ierr); 11918c55e3fSMatthew G. Knepley ierr = VecAXPY(y, -1.0, G);CHKERRQ(ierr); 12018c55e3fSMatthew G. Knepley ierr = TSMimexRestoreXstarAndG(ts, dm, &Xstar, &G);CHKERRQ(ierr); 12118c55e3fSMatthew G. Knepley } 122abadfa0fSMatthew G. Knepley ts->dm = dmsave; 123abadfa0fSMatthew G. Knepley ierr = TSMimexRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 124abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 125abadfa0fSMatthew G. Knepley } 126abadfa0fSMatthew G. Knepley 127abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_Mimex(SNES snes, Vec x, Mat A, Mat B, TS ts) 128abadfa0fSMatthew G. Knepley { 129abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 130abadfa0fSMatthew G. Knepley DM dm, dmsave; 131abadfa0fSMatthew G. Knepley Vec Xdot; 132abadfa0fSMatthew G. Knepley PetscReal shift = 1./ts->time_step; 133abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 134abadfa0fSMatthew G. Knepley 135abadfa0fSMatthew G. Knepley PetscFunctionBegin; 136abadfa0fSMatthew G. Knepley /* th->Xdot has already been computed in SNESTSFormFunction_Mimex (SNES guarantees this) */ 137abadfa0fSMatthew G. Knepley ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 138abadfa0fSMatthew G. Knepley ierr = TSMimexGetX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr); 139abadfa0fSMatthew G. Knepley 140abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 141abadfa0fSMatthew G. Knepley dmsave = ts->dm; 142abadfa0fSMatthew G. Knepley ts->dm = dm; 143abadfa0fSMatthew G. Knepley ierr = TSComputeIJacobian(ts, mimex->stage_time, x, Xdot, shift, A, B, PETSC_TRUE);CHKERRQ(ierr); 144abadfa0fSMatthew G. Knepley ts->dm = dmsave; 145abadfa0fSMatthew G. Knepley ierr = TSMimexRestoreX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr); 146abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 147abadfa0fSMatthew G. Knepley } 148abadfa0fSMatthew G. Knepley 14918c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Split(TS ts) 150abadfa0fSMatthew G. Knepley { 151abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 152abadfa0fSMatthew G. Knepley DM dm; 153abadfa0fSMatthew G. Knepley PetscDS prob; 154abadfa0fSMatthew G. Knepley PetscSection s; 155abadfa0fSMatthew G. Knepley Vec sol = ts->vec_sol, update = mimex->update; 156abadfa0fSMatthew G. Knepley const PetscScalar *aupdate; 157abadfa0fSMatthew G. Knepley PetscScalar *asol, dt = ts->time_step; 158abadfa0fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 159abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 160abadfa0fSMatthew G. Knepley 161abadfa0fSMatthew G. Knepley PetscFunctionBegin; 162abadfa0fSMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 163abadfa0fSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 164*92fd8e1eSJed Brown ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 165abadfa0fSMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 166abadfa0fSMatthew G. Knepley ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 167abadfa0fSMatthew G. Knepley ierr = TSPreStage(ts, ts->ptime);CHKERRQ(ierr); 168abadfa0fSMatthew G. Knepley /* Compute implicit update */ 169abadfa0fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 170abadfa0fSMatthew G. Knepley ierr = VecCopy(sol, update);CHKERRQ(ierr); 171abadfa0fSMatthew G. Knepley ierr = SNESSolve(ts->snes, NULL, update);CHKERRQ(ierr); 172abadfa0fSMatthew G. Knepley ierr = VecGetArrayRead(update, &aupdate);CHKERRQ(ierr); 173abadfa0fSMatthew G. Knepley ierr = VecGetArray(sol, &asol);CHKERRQ(ierr); 174abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 17518c55e3fSMatthew G. Knepley PetscBool implicit; 176abadfa0fSMatthew G. Knepley 17718c55e3fSMatthew G. Knepley ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr); 17818c55e3fSMatthew G. Knepley if (!implicit) continue; 179abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 180abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 181abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 182abadfa0fSMatthew G. Knepley 183abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 184abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr); 185abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au);CHKERRQ(ierr); 186abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRef(dm, p, f, asol, &as);CHKERRQ(ierr); 187abadfa0fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) as[d] = au[d]; 188abadfa0fSMatthew G. Knepley } 189abadfa0fSMatthew G. Knepley } 190abadfa0fSMatthew G. Knepley ierr = VecRestoreArrayRead(update, &aupdate);CHKERRQ(ierr); 191abadfa0fSMatthew G. Knepley ierr = VecRestoreArray(sol, &asol);CHKERRQ(ierr); 192abadfa0fSMatthew G. Knepley /* Compute explicit update */ 193abadfa0fSMatthew G. Knepley ierr = TSComputeRHSFunction(ts, ts->ptime, sol, update);CHKERRQ(ierr); 194abadfa0fSMatthew G. Knepley ierr = VecGetArrayRead(update, &aupdate);CHKERRQ(ierr); 195abadfa0fSMatthew G. Knepley ierr = VecGetArray(sol, &asol);CHKERRQ(ierr); 196abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 19718c55e3fSMatthew G. Knepley PetscBool implicit; 198abadfa0fSMatthew G. Knepley 19918c55e3fSMatthew G. Knepley ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr); 20018c55e3fSMatthew G. Knepley if (implicit) continue; 201abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 202abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 203abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 204abadfa0fSMatthew G. Knepley 205abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 206abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr); 207abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au);CHKERRQ(ierr); 208abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRef(dm, p, f, asol, &as);CHKERRQ(ierr); 209abadfa0fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) as[d] += dt*au[d]; 210abadfa0fSMatthew G. Knepley } 211abadfa0fSMatthew G. Knepley } 212abadfa0fSMatthew G. Knepley ierr = VecRestoreArrayRead(update, &aupdate);CHKERRQ(ierr); 213abadfa0fSMatthew G. Knepley ierr = VecRestoreArray(sol, &asol);CHKERRQ(ierr); 214abadfa0fSMatthew G. Knepley ierr = TSPostStage(ts, ts->ptime, 0, &sol);CHKERRQ(ierr); 215abadfa0fSMatthew G. Knepley ts->ptime += ts->time_step; 216abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 217abadfa0fSMatthew G. Knepley } 218abadfa0fSMatthew G. Knepley 21918c55e3fSMatthew G. Knepley 22018c55e3fSMatthew G. Knepley /* Evalute F at U and G at U0 for explicit fields and U for implicit fields */ 22118c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Implicit(TS ts) 22218c55e3fSMatthew G. Knepley { 22318c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 22418c55e3fSMatthew G. Knepley Vec sol = ts->vec_sol; 22518c55e3fSMatthew G. Knepley Vec update = mimex->update; 22618c55e3fSMatthew G. Knepley PetscErrorCode ierr; 22718c55e3fSMatthew G. Knepley 22818c55e3fSMatthew G. Knepley PetscFunctionBegin; 22918c55e3fSMatthew G. Knepley ierr = TSPreStage(ts, ts->ptime);CHKERRQ(ierr); 23018c55e3fSMatthew G. Knepley /* Compute implicit update */ 23118c55e3fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 23218c55e3fSMatthew G. Knepley ts->ptime += ts->time_step; 23318c55e3fSMatthew G. Knepley ierr = VecCopy(sol, update);CHKERRQ(ierr); 23418c55e3fSMatthew G. Knepley ierr = SNESSolve(ts->snes, NULL, update);CHKERRQ(ierr); 23518c55e3fSMatthew G. Knepley ierr = VecCopy(update, sol);CHKERRQ(ierr); 23618c55e3fSMatthew G. Knepley ierr = TSPostStage(ts, ts->ptime, 0, &sol);CHKERRQ(ierr); 23718c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 23818c55e3fSMatthew G. Knepley } 23918c55e3fSMatthew G. Knepley 24018c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex(TS ts) 24118c55e3fSMatthew G. Knepley { 24218c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 24318c55e3fSMatthew G. Knepley PetscErrorCode ierr; 24418c55e3fSMatthew G. Knepley 24518c55e3fSMatthew G. Knepley PetscFunctionBegin; 24618c55e3fSMatthew G. Knepley switch(mimex->version) { 24718c55e3fSMatthew G. Knepley case 0: 24818c55e3fSMatthew G. Knepley ierr = TSStep_Mimex_Split(ts);CHKERRQ(ierr); break; 24918c55e3fSMatthew G. Knepley case 1: 25018c55e3fSMatthew G. Knepley ierr = TSStep_Mimex_Implicit(ts);CHKERRQ(ierr); break; 25118c55e3fSMatthew G. Knepley default: 25218c55e3fSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_OUTOFRANGE, "Unknown MIMEX version %d", mimex->version); 25318c55e3fSMatthew G. Knepley } 25418c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 25518c55e3fSMatthew G. Knepley } 25618c55e3fSMatthew G. Knepley 257abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 258abadfa0fSMatthew G. Knepley 259abadfa0fSMatthew G. Knepley static PetscErrorCode TSSetUp_Mimex(TS ts) 260abadfa0fSMatthew G. Knepley { 261abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 262abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 263abadfa0fSMatthew G. Knepley 264abadfa0fSMatthew G. Knepley PetscFunctionBegin; 265abadfa0fSMatthew G. Knepley ierr = VecDuplicate(ts->vec_sol, &mimex->update);CHKERRQ(ierr); 266abadfa0fSMatthew G. Knepley ierr = VecDuplicate(ts->vec_sol, &mimex->Xdot);CHKERRQ(ierr); 267abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 268abadfa0fSMatthew G. Knepley } 269abadfa0fSMatthew G. Knepley 270abadfa0fSMatthew G. Knepley static PetscErrorCode TSReset_Mimex(TS ts) 271abadfa0fSMatthew G. Knepley { 272abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 273abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 274abadfa0fSMatthew G. Knepley 275abadfa0fSMatthew G. Knepley PetscFunctionBegin; 276abadfa0fSMatthew G. Knepley ierr = VecDestroy(&mimex->update);CHKERRQ(ierr); 277abadfa0fSMatthew G. Knepley ierr = VecDestroy(&mimex->Xdot);CHKERRQ(ierr); 278abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 279abadfa0fSMatthew G. Knepley } 280abadfa0fSMatthew G. Knepley 281abadfa0fSMatthew G. Knepley static PetscErrorCode TSDestroy_Mimex(TS ts) 282abadfa0fSMatthew G. Knepley { 283abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 284abadfa0fSMatthew G. Knepley 285abadfa0fSMatthew G. Knepley PetscFunctionBegin; 286abadfa0fSMatthew G. Knepley ierr = TSReset_Mimex(ts);CHKERRQ(ierr); 287abadfa0fSMatthew G. Knepley ierr = PetscFree(ts->data);CHKERRQ(ierr); 288abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 289abadfa0fSMatthew G. Knepley } 290abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 291abadfa0fSMatthew G. Knepley 2924416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Mimex(PetscOptionItems *PetscOptionsObject, TS ts) 293abadfa0fSMatthew G. Knepley { 29418c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 29518c55e3fSMatthew G. Knepley PetscErrorCode ierr; 29618c55e3fSMatthew G. Knepley 297abadfa0fSMatthew G. Knepley PetscFunctionBegin; 29818c55e3fSMatthew G. Knepley ierr = PetscOptionsHead(PetscOptionsObject, "MIMEX ODE solver options");CHKERRQ(ierr); 29918c55e3fSMatthew G. Knepley { 30018c55e3fSMatthew G. Knepley ierr = PetscOptionsInt("-ts_mimex_version", "Algorithm version", "TSMimexSetVersion", mimex->version, &mimex->version, NULL);CHKERRQ(ierr); 30118c55e3fSMatthew G. Knepley } 30218c55e3fSMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 303abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 304abadfa0fSMatthew G. Knepley } 305abadfa0fSMatthew G. Knepley 306abadfa0fSMatthew G. Knepley static PetscErrorCode TSView_Mimex(TS ts,PetscViewer viewer) 307abadfa0fSMatthew G. Knepley { 3085bf398f3SMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 3095bf398f3SMatthew G. Knepley PetscBool iascii; 3105bf398f3SMatthew G. Knepley PetscErrorCode ierr; 3115bf398f3SMatthew G. Knepley 312abadfa0fSMatthew G. Knepley PetscFunctionBegin; 3135bf398f3SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 3145bf398f3SMatthew G. Knepley if (iascii) { 3155bf398f3SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, " Version = %D\n", mimex->version);CHKERRQ(ierr); 3165bf398f3SMatthew G. Knepley } 317abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 318abadfa0fSMatthew G. Knepley } 319abadfa0fSMatthew G. Knepley 320abadfa0fSMatthew G. Knepley static PetscErrorCode TSInterpolate_Mimex(TS ts,PetscReal t,Vec X) 321abadfa0fSMatthew G. Knepley { 322abadfa0fSMatthew G. Knepley PetscReal alpha = (ts->ptime - t)/ts->time_step; 323abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 324abadfa0fSMatthew G. Knepley 325abadfa0fSMatthew G. Knepley PetscFunctionBegin; 326abadfa0fSMatthew G. Knepley ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr); 327abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 328abadfa0fSMatthew G. Knepley } 329abadfa0fSMatthew G. Knepley 330560360afSLisandro Dalcin static PetscErrorCode TSComputeLinearStability_Mimex(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 331abadfa0fSMatthew G. Knepley { 332abadfa0fSMatthew G. Knepley PetscFunctionBegin; 333abadfa0fSMatthew G. Knepley *yr = 1.0 + xr; 334abadfa0fSMatthew G. Knepley *yi = xi; 335abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 336abadfa0fSMatthew G. Knepley } 337abadfa0fSMatthew G. Knepley /* ------------------------------------------------------------ */ 338abadfa0fSMatthew G. Knepley 339abadfa0fSMatthew G. Knepley /*MC 340abadfa0fSMatthew G. Knepley TSMIMEX - ODE solver using the explicit forward Mimex method 341abadfa0fSMatthew G. Knepley 342abadfa0fSMatthew G. Knepley Level: beginner 343abadfa0fSMatthew G. Knepley 344abadfa0fSMatthew G. Knepley .seealso: TSCreate(), TS, TSSetType(), TSBEULER 345abadfa0fSMatthew G. Knepley 346abadfa0fSMatthew G. Knepley M*/ 347abadfa0fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_Mimex(TS ts) 348abadfa0fSMatthew G. Knepley { 349abadfa0fSMatthew G. Knepley TS_Mimex *mimex; 350abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 351abadfa0fSMatthew G. Knepley 352abadfa0fSMatthew G. Knepley PetscFunctionBegin; 353abadfa0fSMatthew G. Knepley ts->ops->setup = TSSetUp_Mimex; 354abadfa0fSMatthew G. Knepley ts->ops->step = TSStep_Mimex; 355abadfa0fSMatthew G. Knepley ts->ops->reset = TSReset_Mimex; 356abadfa0fSMatthew G. Knepley ts->ops->destroy = TSDestroy_Mimex; 357abadfa0fSMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_Mimex; 358abadfa0fSMatthew G. Knepley ts->ops->view = TSView_Mimex; 359abadfa0fSMatthew G. Knepley ts->ops->interpolate = TSInterpolate_Mimex; 360abadfa0fSMatthew G. Knepley ts->ops->linearstability = TSComputeLinearStability_Mimex; 361abadfa0fSMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_Mimex; 362abadfa0fSMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_Mimex; 3632ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 364abadfa0fSMatthew G. Knepley 365abadfa0fSMatthew G. Knepley ierr = PetscNewLog(ts,&mimex);CHKERRQ(ierr); 366abadfa0fSMatthew G. Knepley ts->data = (void*)mimex; 367fdc0a3ceSMatthew G. Knepley 368fdc0a3ceSMatthew G. Knepley mimex->version = 1; 369abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 370abadfa0fSMatthew G. Knepley } 371