1abadfa0fSMatthew G. Knepley /* 2abadfa0fSMatthew G. Knepley Code for Timestepping with my makeshift IMEX. 3abadfa0fSMatthew G. Knepley */ 4abadfa0fSMatthew G. Knepley #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 { 9*18c55e3fSMatthew G. Knepley Vec Xdot, update, Xstar, G; 10abadfa0fSMatthew G. Knepley PetscReal stage_time; 11*18c55e3fSMatthew G. Knepley PetscInt version; 12abadfa0fSMatthew G. Knepley } TS_Mimex; 13abadfa0fSMatthew G. Knepley 14abadfa0fSMatthew G. Knepley #undef __FUNCT__ 15abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSMimexGetX0AndXdot" 16abadfa0fSMatthew G. Knepley static PetscErrorCode TSMimexGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 17abadfa0fSMatthew G. Knepley { 18abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 19abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 20abadfa0fSMatthew G. Knepley 21abadfa0fSMatthew G. Knepley PetscFunctionBegin; 22abadfa0fSMatthew G. Knepley if (X0) { 23abadfa0fSMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSMimex_X0", X0);CHKERRQ(ierr);} 24abadfa0fSMatthew G. Knepley else {*X0 = ts->vec_sol;} 25abadfa0fSMatthew G. Knepley } 26abadfa0fSMatthew G. Knepley if (Xdot) { 27abadfa0fSMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSMimex_Xdot", Xdot);CHKERRQ(ierr);} 28abadfa0fSMatthew G. Knepley else {*Xdot = mimex->Xdot;} 29abadfa0fSMatthew G. Knepley } 30abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 31abadfa0fSMatthew G. Knepley } 32abadfa0fSMatthew G. Knepley 33abadfa0fSMatthew G. Knepley #undef __FUNCT__ 34abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSMimexRestoreX0AndXdot" 35abadfa0fSMatthew G. Knepley static PetscErrorCode TSMimexRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 36abadfa0fSMatthew G. Knepley { 37abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 38abadfa0fSMatthew G. Knepley 39abadfa0fSMatthew G. Knepley PetscFunctionBegin; 40abadfa0fSMatthew G. Knepley if (X0) if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_X0", X0);CHKERRQ(ierr);} 41abadfa0fSMatthew G. Knepley if (Xdot) if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_Xdot", Xdot);CHKERRQ(ierr);} 42abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 43abadfa0fSMatthew G. Knepley } 44abadfa0fSMatthew G. Knepley 45*18c55e3fSMatthew G. Knepley #undef __FUNCT__ 46*18c55e3fSMatthew G. Knepley #define __FUNCT__ "TSMimexGetXstarAndG" 47*18c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexGetXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 48*18c55e3fSMatthew G. Knepley { 49*18c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 50*18c55e3fSMatthew G. Knepley PetscErrorCode ierr; 51*18c55e3fSMatthew G. Knepley 52*18c55e3fSMatthew G. Knepley PetscFunctionBegin; 53*18c55e3fSMatthew G. Knepley if (!mimex->Xstar) {ierr = VecDuplicate(ts->vec_sol, &mimex->Xstar);CHKERRQ(ierr);} 54*18c55e3fSMatthew G. Knepley if (!mimex->G) {ierr = VecDuplicate(ts->vec_sol, &mimex->G);CHKERRQ(ierr);} 55*18c55e3fSMatthew G. Knepley *Xstar = mimex->Xstar; 56*18c55e3fSMatthew G. Knepley *G = mimex->G; 57*18c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 58*18c55e3fSMatthew G. Knepley } 59*18c55e3fSMatthew G. Knepley 60*18c55e3fSMatthew G. Knepley #undef __FUNCT__ 61*18c55e3fSMatthew G. Knepley #define __FUNCT__ "TSMimexRestoreXstarAndG" 62*18c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexRestoreXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G) 63*18c55e3fSMatthew G. Knepley { 64*18c55e3fSMatthew G. Knepley PetscFunctionBegin; 65*18c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 66*18c55e3fSMatthew G. Knepley } 67*18c55e3fSMatthew G. Knepley 68abadfa0fSMatthew G. Knepley /* 69abadfa0fSMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 70abadfa0fSMatthew G. Knepley G(U) = F[t0+dt, U, (U-U0)*shift] = 0 71abadfa0fSMatthew G. Knepley */ 72abadfa0fSMatthew G. Knepley #undef __FUNCT__ 73abadfa0fSMatthew G. Knepley #define __FUNCT__ "SNESTSFormFunction_Mimex" 74abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_Mimex(SNES snes, Vec x, Vec y, TS ts) 75abadfa0fSMatthew G. Knepley { 76abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 77abadfa0fSMatthew G. Knepley DM dm, dmsave; 78abadfa0fSMatthew G. Knepley Vec X0, Xdot; 79abadfa0fSMatthew G. Knepley PetscReal shift = 1./ts->time_step; 80abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 81abadfa0fSMatthew G. Knepley 82abadfa0fSMatthew G. Knepley PetscFunctionBegin; 83abadfa0fSMatthew G. Knepley ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 84abadfa0fSMatthew G. Knepley ierr = TSMimexGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 85abadfa0fSMatthew G. Knepley ierr = VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x);CHKERRQ(ierr); 86abadfa0fSMatthew G. Knepley 87abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 88abadfa0fSMatthew G. Knepley dmsave = ts->dm; 89abadfa0fSMatthew G. Knepley ts->dm = dm; 90abadfa0fSMatthew G. Knepley ierr = TSComputeIFunction(ts, mimex->stage_time, x, Xdot, y, PETSC_TRUE);CHKERRQ(ierr); 91*18c55e3fSMatthew G. Knepley if (mimex->version == 1) { 92*18c55e3fSMatthew G. Knepley DM dm; 93*18c55e3fSMatthew G. Knepley PetscDS prob; 94*18c55e3fSMatthew G. Knepley PetscSection s; 95*18c55e3fSMatthew G. Knepley Vec Xstar, G; 96*18c55e3fSMatthew G. Knepley const PetscScalar *ax; 97*18c55e3fSMatthew G. Knepley PetscScalar *axstar; 98*18c55e3fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 99*18c55e3fSMatthew G. Knepley 100*18c55e3fSMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 101*18c55e3fSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 102*18c55e3fSMatthew G. Knepley ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 103*18c55e3fSMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 104*18c55e3fSMatthew G. Knepley ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 105*18c55e3fSMatthew G. Knepley ierr = TSMimexGetXstarAndG(ts, dm, &Xstar, &G);CHKERRQ(ierr); 106*18c55e3fSMatthew G. Knepley ierr = VecCopy(X0, Xstar);CHKERRQ(ierr); 107*18c55e3fSMatthew G. Knepley ierr = VecGetArrayRead(x, &ax);CHKERRQ(ierr); 108*18c55e3fSMatthew G. Knepley ierr = VecGetArray(Xstar, &axstar);CHKERRQ(ierr); 109*18c55e3fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 110*18c55e3fSMatthew G. Knepley PetscBool implicit; 111*18c55e3fSMatthew G. Knepley 112*18c55e3fSMatthew G. Knepley ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr); 113*18c55e3fSMatthew G. Knepley if (!implicit) continue; 114*18c55e3fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 115*18c55e3fSMatthew G. Knepley PetscScalar *a, *axs; 116*18c55e3fSMatthew G. Knepley PetscInt fdof, fcdof, d; 117*18c55e3fSMatthew G. Knepley 118*18c55e3fSMatthew G. Knepley ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 119*18c55e3fSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr); 120*18c55e3fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRead(dm, p, f, ax, &a);CHKERRQ(ierr); 121*18c55e3fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRef(dm, p, f, axstar, &axs);CHKERRQ(ierr); 122*18c55e3fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) axs[d] = a[d]; 123*18c55e3fSMatthew G. Knepley } 124*18c55e3fSMatthew G. Knepley } 125*18c55e3fSMatthew G. Knepley ierr = VecRestoreArrayRead(x, &ax);CHKERRQ(ierr); 126*18c55e3fSMatthew G. Knepley ierr = VecRestoreArray(Xstar, &axstar);CHKERRQ(ierr); 127*18c55e3fSMatthew G. Knepley ierr = TSComputeRHSFunction(ts, ts->ptime, Xstar, G);CHKERRQ(ierr); 128*18c55e3fSMatthew G. Knepley ierr = VecAXPY(y, -1.0, G);CHKERRQ(ierr); 129*18c55e3fSMatthew G. Knepley ierr = TSMimexRestoreXstarAndG(ts, dm, &Xstar, &G);CHKERRQ(ierr); 130*18c55e3fSMatthew G. Knepley } 131abadfa0fSMatthew G. Knepley ts->dm = dmsave; 132abadfa0fSMatthew G. Knepley ierr = TSMimexRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 133abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 134abadfa0fSMatthew G. Knepley } 135abadfa0fSMatthew G. Knepley 136abadfa0fSMatthew G. Knepley #undef __FUNCT__ 137abadfa0fSMatthew G. Knepley #define __FUNCT__ "SNESTSFormJacobian_Mimex" 138abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_Mimex(SNES snes, Vec x, Mat A, Mat B, TS ts) 139abadfa0fSMatthew G. Knepley { 140abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 141abadfa0fSMatthew G. Knepley DM dm, dmsave; 142abadfa0fSMatthew G. Knepley Vec Xdot; 143abadfa0fSMatthew G. Knepley PetscReal shift = 1./ts->time_step; 144abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 145abadfa0fSMatthew G. Knepley 146abadfa0fSMatthew G. Knepley PetscFunctionBegin; 147abadfa0fSMatthew G. Knepley /* th->Xdot has already been computed in SNESTSFormFunction_Mimex (SNES guarantees this) */ 148abadfa0fSMatthew G. Knepley ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 149abadfa0fSMatthew G. Knepley ierr = TSMimexGetX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr); 150abadfa0fSMatthew G. Knepley 151abadfa0fSMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 152abadfa0fSMatthew G. Knepley dmsave = ts->dm; 153abadfa0fSMatthew G. Knepley ts->dm = dm; 154abadfa0fSMatthew G. Knepley ierr = TSComputeIJacobian(ts, mimex->stage_time, x, Xdot, shift, A, B, PETSC_TRUE);CHKERRQ(ierr); 155abadfa0fSMatthew G. Knepley ts->dm = dmsave; 156abadfa0fSMatthew G. Knepley ierr = TSMimexRestoreX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr); 157abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 158abadfa0fSMatthew G. Knepley } 159abadfa0fSMatthew G. Knepley 160abadfa0fSMatthew G. Knepley #undef __FUNCT__ 161*18c55e3fSMatthew G. Knepley #define __FUNCT__ "TSStep_Mimex_Split" 162*18c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Split(TS ts) 163abadfa0fSMatthew G. Knepley { 164abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 165abadfa0fSMatthew G. Knepley DM dm; 166abadfa0fSMatthew G. Knepley PetscDS prob; 167abadfa0fSMatthew G. Knepley PetscSection s; 168abadfa0fSMatthew G. Knepley Vec sol = ts->vec_sol, update = mimex->update; 169abadfa0fSMatthew G. Knepley const PetscScalar *aupdate; 170abadfa0fSMatthew G. Knepley PetscScalar *asol, dt = ts->time_step; 171abadfa0fSMatthew G. Knepley PetscInt Nf, f, pStart, pEnd, p; 172abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 173abadfa0fSMatthew G. Knepley 174abadfa0fSMatthew G. Knepley PetscFunctionBegin; 175abadfa0fSMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 176abadfa0fSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 177abadfa0fSMatthew G. Knepley ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 178abadfa0fSMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 179abadfa0fSMatthew G. Knepley ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 180abadfa0fSMatthew G. Knepley ierr = TSPreStep(ts);CHKERRQ(ierr); 181abadfa0fSMatthew G. Knepley ierr = TSPreStage(ts, ts->ptime);CHKERRQ(ierr); 182abadfa0fSMatthew G. Knepley /* Compute implicit update */ 183abadfa0fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 184abadfa0fSMatthew G. Knepley ierr = VecCopy(sol, update);CHKERRQ(ierr); 185abadfa0fSMatthew G. Knepley ierr = SNESSolve(ts->snes, NULL, update);CHKERRQ(ierr); 186abadfa0fSMatthew G. Knepley ierr = VecGetArrayRead(update, &aupdate);CHKERRQ(ierr); 187abadfa0fSMatthew G. Knepley ierr = VecGetArray(sol, &asol);CHKERRQ(ierr); 188abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 189*18c55e3fSMatthew G. Knepley PetscBool implicit; 190abadfa0fSMatthew G. Knepley 191*18c55e3fSMatthew G. Knepley ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr); 192*18c55e3fSMatthew G. Knepley if (!implicit) continue; 193abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 194abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 195abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 196abadfa0fSMatthew G. Knepley 197abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 198abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr); 199abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au);CHKERRQ(ierr); 200abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRef(dm, p, f, asol, &as);CHKERRQ(ierr); 201abadfa0fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) as[d] = au[d]; 202abadfa0fSMatthew G. Knepley } 203abadfa0fSMatthew G. Knepley } 204abadfa0fSMatthew G. Knepley ierr = VecRestoreArrayRead(update, &aupdate);CHKERRQ(ierr); 205abadfa0fSMatthew G. Knepley ierr = VecRestoreArray(sol, &asol);CHKERRQ(ierr); 206abadfa0fSMatthew G. Knepley /* Compute explicit update */ 207abadfa0fSMatthew G. Knepley ierr = TSComputeRHSFunction(ts, ts->ptime, sol, update);CHKERRQ(ierr); 208abadfa0fSMatthew G. Knepley ierr = VecGetArrayRead(update, &aupdate);CHKERRQ(ierr); 209abadfa0fSMatthew G. Knepley ierr = VecGetArray(sol, &asol);CHKERRQ(ierr); 210abadfa0fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 211*18c55e3fSMatthew G. Knepley PetscBool implicit; 212abadfa0fSMatthew G. Knepley 213*18c55e3fSMatthew G. Knepley ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr); 214*18c55e3fSMatthew G. Knepley if (implicit) continue; 215abadfa0fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 216abadfa0fSMatthew G. Knepley PetscScalar *au, *as; 217abadfa0fSMatthew G. Knepley PetscInt fdof, fcdof, d; 218abadfa0fSMatthew G. Knepley 219abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 220abadfa0fSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr); 221abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au);CHKERRQ(ierr); 222abadfa0fSMatthew G. Knepley ierr = DMPlexPointGlobalFieldRef(dm, p, f, asol, &as);CHKERRQ(ierr); 223abadfa0fSMatthew G. Knepley for (d = 0; d < fdof-fcdof; ++d) as[d] += dt*au[d]; 224abadfa0fSMatthew G. Knepley } 225abadfa0fSMatthew G. Knepley } 226abadfa0fSMatthew G. Knepley ierr = VecRestoreArrayRead(update, &aupdate);CHKERRQ(ierr); 227abadfa0fSMatthew G. Knepley ierr = VecRestoreArray(sol, &asol);CHKERRQ(ierr); 228abadfa0fSMatthew G. Knepley ierr = TSPostStage(ts, ts->ptime, 0, &sol);CHKERRQ(ierr); 229abadfa0fSMatthew G. Knepley ts->ptime += ts->time_step; 230abadfa0fSMatthew G. Knepley ts->steps++; 231abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 232abadfa0fSMatthew G. Knepley } 233abadfa0fSMatthew G. Knepley 234*18c55e3fSMatthew G. Knepley 235*18c55e3fSMatthew G. Knepley #undef __FUNCT__ 236*18c55e3fSMatthew G. Knepley #define __FUNCT__ "TSStep_Mimex_Implicit" 237*18c55e3fSMatthew G. Knepley /* Evalute F at U and G at U0 for explicit fields and U for implicit fields */ 238*18c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Implicit(TS ts) 239*18c55e3fSMatthew G. Knepley { 240*18c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 241*18c55e3fSMatthew G. Knepley Vec sol = ts->vec_sol; 242*18c55e3fSMatthew G. Knepley Vec update = mimex->update; 243*18c55e3fSMatthew G. Knepley PetscErrorCode ierr; 244*18c55e3fSMatthew G. Knepley 245*18c55e3fSMatthew G. Knepley PetscFunctionBegin; 246*18c55e3fSMatthew G. Knepley ierr = TSPreStep(ts);CHKERRQ(ierr); 247*18c55e3fSMatthew G. Knepley ierr = TSPreStage(ts, ts->ptime);CHKERRQ(ierr); 248*18c55e3fSMatthew G. Knepley /* Compute implicit update */ 249*18c55e3fSMatthew G. Knepley mimex->stage_time = ts->ptime + ts->time_step; 250*18c55e3fSMatthew G. Knepley ts->ptime += ts->time_step; 251*18c55e3fSMatthew G. Knepley ierr = VecCopy(sol, update);CHKERRQ(ierr); 252*18c55e3fSMatthew G. Knepley ierr = SNESSolve(ts->snes, NULL, update);CHKERRQ(ierr); 253*18c55e3fSMatthew G. Knepley ierr = VecCopy(update, sol);CHKERRQ(ierr); 254*18c55e3fSMatthew G. Knepley ierr = TSPostStage(ts, ts->ptime, 0, &sol);CHKERRQ(ierr); 255*18c55e3fSMatthew G. Knepley ts->steps++; 256*18c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 257*18c55e3fSMatthew G. Knepley } 258*18c55e3fSMatthew G. Knepley 259*18c55e3fSMatthew G. Knepley #undef __FUNCT__ 260*18c55e3fSMatthew G. Knepley #define __FUNCT__ "TSStep_Mimex" 261*18c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex(TS ts) 262*18c55e3fSMatthew G. Knepley { 263*18c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 264*18c55e3fSMatthew G. Knepley PetscErrorCode ierr; 265*18c55e3fSMatthew G. Knepley 266*18c55e3fSMatthew G. Knepley PetscFunctionBegin; 267*18c55e3fSMatthew G. Knepley switch(mimex->version) { 268*18c55e3fSMatthew G. Knepley case 0: 269*18c55e3fSMatthew G. Knepley ierr = TSStep_Mimex_Split(ts);CHKERRQ(ierr); break; 270*18c55e3fSMatthew G. Knepley case 1: 271*18c55e3fSMatthew G. Knepley ierr = TSStep_Mimex_Implicit(ts);CHKERRQ(ierr); break; 272*18c55e3fSMatthew G. Knepley default: 273*18c55e3fSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_OUTOFRANGE, "Unknown MIMEX version %d", mimex->version); 274*18c55e3fSMatthew G. Knepley } 275*18c55e3fSMatthew G. Knepley PetscFunctionReturn(0); 276*18c55e3fSMatthew G. Knepley } 277*18c55e3fSMatthew G. Knepley 278abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 279abadfa0fSMatthew G. Knepley 280abadfa0fSMatthew G. Knepley #undef __FUNCT__ 281abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSSetUp_Mimex" 282abadfa0fSMatthew G. Knepley static PetscErrorCode TSSetUp_Mimex(TS ts) 283abadfa0fSMatthew G. Knepley { 284abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 285abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 286abadfa0fSMatthew G. Knepley 287abadfa0fSMatthew G. Knepley PetscFunctionBegin; 288abadfa0fSMatthew G. Knepley ierr = VecDuplicate(ts->vec_sol, &mimex->update);CHKERRQ(ierr); 289abadfa0fSMatthew G. Knepley ierr = VecDuplicate(ts->vec_sol, &mimex->Xdot);CHKERRQ(ierr); 290abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 291abadfa0fSMatthew G. Knepley } 292abadfa0fSMatthew G. Knepley 293abadfa0fSMatthew G. Knepley #undef __FUNCT__ 294abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSReset_Mimex" 295abadfa0fSMatthew G. Knepley static PetscErrorCode TSReset_Mimex(TS ts) 296abadfa0fSMatthew G. Knepley { 297abadfa0fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex*)ts->data; 298abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 299abadfa0fSMatthew G. Knepley 300abadfa0fSMatthew G. Knepley PetscFunctionBegin; 301*18c55e3fSMatthew G. Knepley ierr = VecDestroy(&mimex->Xstar);CHKERRQ(ierr); 302*18c55e3fSMatthew G. Knepley ierr = VecDestroy(&mimex->G);CHKERRQ(ierr); 303abadfa0fSMatthew G. Knepley ierr = VecDestroy(&mimex->update);CHKERRQ(ierr); 304abadfa0fSMatthew G. Knepley ierr = VecDestroy(&mimex->Xdot);CHKERRQ(ierr); 305abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 306abadfa0fSMatthew G. Knepley } 307abadfa0fSMatthew G. Knepley 308abadfa0fSMatthew G. Knepley #undef __FUNCT__ 309abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSDestroy_Mimex" 310abadfa0fSMatthew G. Knepley static PetscErrorCode TSDestroy_Mimex(TS ts) 311abadfa0fSMatthew G. Knepley { 312abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 313abadfa0fSMatthew G. Knepley 314abadfa0fSMatthew G. Knepley PetscFunctionBegin; 315abadfa0fSMatthew G. Knepley ierr = TSReset_Mimex(ts);CHKERRQ(ierr); 316abadfa0fSMatthew G. Knepley ierr = PetscFree(ts->data);CHKERRQ(ierr); 317abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 318abadfa0fSMatthew G. Knepley } 319abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/ 320abadfa0fSMatthew G. Knepley 321abadfa0fSMatthew G. Knepley #undef __FUNCT__ 322abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSSetFromOptions_Mimex" 323abadfa0fSMatthew G. Knepley static PetscErrorCode TSSetFromOptions_Mimex(PetscOptions *PetscOptionsObject, TS ts) 324abadfa0fSMatthew G. Knepley { 325*18c55e3fSMatthew G. Knepley TS_Mimex *mimex = (TS_Mimex *) ts->data; 326*18c55e3fSMatthew G. Knepley PetscErrorCode ierr; 327*18c55e3fSMatthew G. Knepley 328abadfa0fSMatthew G. Knepley PetscFunctionBegin; 329*18c55e3fSMatthew G. Knepley ierr = PetscOptionsHead(PetscOptionsObject, "MIMEX ODE solver options");CHKERRQ(ierr); 330*18c55e3fSMatthew G. Knepley { 331*18c55e3fSMatthew G. Knepley ierr = PetscOptionsInt("-ts_mimex_version", "Algorithm version", "TSMimexSetVersion", mimex->version, &mimex->version, NULL);CHKERRQ(ierr); 332*18c55e3fSMatthew G. Knepley } 333*18c55e3fSMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 334abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 335abadfa0fSMatthew G. Knepley } 336abadfa0fSMatthew G. Knepley 337abadfa0fSMatthew G. Knepley #undef __FUNCT__ 338abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSView_Mimex" 339abadfa0fSMatthew G. Knepley static PetscErrorCode TSView_Mimex(TS ts,PetscViewer viewer) 340abadfa0fSMatthew G. Knepley { 341abadfa0fSMatthew G. Knepley PetscFunctionBegin; 342abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 343abadfa0fSMatthew G. Knepley } 344abadfa0fSMatthew G. Knepley 345abadfa0fSMatthew G. Knepley #undef __FUNCT__ 346abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSInterpolate_Mimex" 347abadfa0fSMatthew G. Knepley static PetscErrorCode TSInterpolate_Mimex(TS ts,PetscReal t,Vec X) 348abadfa0fSMatthew G. Knepley { 349abadfa0fSMatthew G. Knepley PetscReal alpha = (ts->ptime - t)/ts->time_step; 350abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 351abadfa0fSMatthew G. Knepley 352abadfa0fSMatthew G. Knepley PetscFunctionBegin; 353abadfa0fSMatthew G. Knepley ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr); 354abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 355abadfa0fSMatthew G. Knepley } 356abadfa0fSMatthew G. Knepley 357abadfa0fSMatthew G. Knepley #undef __FUNCT__ 358abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSComputeLinearStability_Mimex" 359abadfa0fSMatthew G. Knepley PetscErrorCode TSComputeLinearStability_Mimex(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 360abadfa0fSMatthew G. Knepley { 361abadfa0fSMatthew G. Knepley PetscFunctionBegin; 362abadfa0fSMatthew G. Knepley *yr = 1.0 + xr; 363abadfa0fSMatthew G. Knepley *yi = xi; 364abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 365abadfa0fSMatthew G. Knepley } 366abadfa0fSMatthew G. Knepley /* ------------------------------------------------------------ */ 367abadfa0fSMatthew G. Knepley 368abadfa0fSMatthew G. Knepley /*MC 369abadfa0fSMatthew G. Knepley TSMIMEX - ODE solver using the explicit forward Mimex method 370abadfa0fSMatthew G. Knepley 371abadfa0fSMatthew G. Knepley Level: beginner 372abadfa0fSMatthew G. Knepley 373abadfa0fSMatthew G. Knepley .seealso: TSCreate(), TS, TSSetType(), TSBEULER 374abadfa0fSMatthew G. Knepley 375abadfa0fSMatthew G. Knepley M*/ 376abadfa0fSMatthew G. Knepley #undef __FUNCT__ 377abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSCreate_Mimex" 378abadfa0fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_Mimex(TS ts) 379abadfa0fSMatthew G. Knepley { 380abadfa0fSMatthew G. Knepley TS_Mimex *mimex; 381abadfa0fSMatthew G. Knepley PetscErrorCode ierr; 382abadfa0fSMatthew G. Knepley 383abadfa0fSMatthew G. Knepley PetscFunctionBegin; 384abadfa0fSMatthew G. Knepley ts->ops->setup = TSSetUp_Mimex; 385abadfa0fSMatthew G. Knepley ts->ops->step = TSStep_Mimex; 386abadfa0fSMatthew G. Knepley ts->ops->reset = TSReset_Mimex; 387abadfa0fSMatthew G. Knepley ts->ops->destroy = TSDestroy_Mimex; 388abadfa0fSMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_Mimex; 389abadfa0fSMatthew G. Knepley ts->ops->view = TSView_Mimex; 390abadfa0fSMatthew G. Knepley ts->ops->interpolate = TSInterpolate_Mimex; 391abadfa0fSMatthew G. Knepley ts->ops->linearstability = TSComputeLinearStability_Mimex; 392abadfa0fSMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_Mimex; 393abadfa0fSMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_Mimex; 394abadfa0fSMatthew G. Knepley 395abadfa0fSMatthew G. Knepley ierr = PetscNewLog(ts,&mimex);CHKERRQ(ierr); 396abadfa0fSMatthew G. Knepley ts->data = (void*)mimex; 397abadfa0fSMatthew G. Knepley PetscFunctionReturn(0); 398abadfa0fSMatthew G. Knepley } 399