xref: /petsc/src/ts/impls/mimex/mimex.c (revision 560360af1e0197128c3ad6271bbfa2e76ad345d4)
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 #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 
4518c55e3fSMatthew G. Knepley #undef __FUNCT__
4618c55e3fSMatthew G. Knepley #define __FUNCT__ "TSMimexGetXstarAndG"
4718c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexGetXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G)
4818c55e3fSMatthew G. Knepley {
4918c55e3fSMatthew G. Knepley   PetscErrorCode ierr;
5018c55e3fSMatthew G. Knepley 
5118c55e3fSMatthew G. Knepley   PetscFunctionBegin;
52010837a9SMatthew G. Knepley   ierr = DMGetNamedGlobalVector(dm, "TSMimex_Xstar", Xstar);CHKERRQ(ierr);
53010837a9SMatthew G. Knepley   ierr = DMGetNamedGlobalVector(dm, "TSMimex_G", G);CHKERRQ(ierr);
5418c55e3fSMatthew G. Knepley   PetscFunctionReturn(0);
5518c55e3fSMatthew G. Knepley }
5618c55e3fSMatthew G. Knepley 
5718c55e3fSMatthew G. Knepley #undef __FUNCT__
5818c55e3fSMatthew G. Knepley #define __FUNCT__ "TSMimexRestoreXstarAndG"
5918c55e3fSMatthew G. Knepley static PetscErrorCode TSMimexRestoreXstarAndG(TS ts, DM dm, Vec *Xstar, Vec *G)
6018c55e3fSMatthew G. Knepley {
61010837a9SMatthew G. Knepley   PetscErrorCode ierr;
62010837a9SMatthew G. Knepley 
6318c55e3fSMatthew G. Knepley   PetscFunctionBegin;
64010837a9SMatthew G. Knepley   ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_Xstar", Xstar);CHKERRQ(ierr);
65010837a9SMatthew G. Knepley   ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_G", G);CHKERRQ(ierr);
6618c55e3fSMatthew G. Knepley   PetscFunctionReturn(0);
6718c55e3fSMatthew G. Knepley }
6818c55e3fSMatthew G. Knepley 
69abadfa0fSMatthew G. Knepley /*
70abadfa0fSMatthew G. Knepley   This defines the nonlinear equation that is to be solved with SNES
71abadfa0fSMatthew G. Knepley   G(U) = F[t0+dt, U, (U-U0)*shift] = 0
72abadfa0fSMatthew G. Knepley */
73abadfa0fSMatthew G. Knepley #undef __FUNCT__
74abadfa0fSMatthew G. Knepley #define __FUNCT__ "SNESTSFormFunction_Mimex"
75abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_Mimex(SNES snes, Vec x, Vec y, TS ts)
76abadfa0fSMatthew G. Knepley {
77abadfa0fSMatthew G. Knepley   TS_Mimex      *mimex = (TS_Mimex *) ts->data;
78abadfa0fSMatthew G. Knepley   DM             dm, dmsave;
79abadfa0fSMatthew G. Knepley   Vec            X0, Xdot;
80abadfa0fSMatthew G. Knepley   PetscReal      shift = 1./ts->time_step;
81abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
82abadfa0fSMatthew G. Knepley 
83abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
84abadfa0fSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
85abadfa0fSMatthew G. Knepley   ierr = TSMimexGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
86abadfa0fSMatthew G. Knepley   ierr = VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x);CHKERRQ(ierr);
87abadfa0fSMatthew G. Knepley 
88abadfa0fSMatthew G. Knepley   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
89abadfa0fSMatthew G. Knepley   dmsave = ts->dm;
90abadfa0fSMatthew G. Knepley   ts->dm = dm;
91abadfa0fSMatthew G. Knepley   ierr   = TSComputeIFunction(ts, mimex->stage_time, x, Xdot, y, PETSC_TRUE);CHKERRQ(ierr);
9218c55e3fSMatthew G. Knepley   if (mimex->version == 1) {
9318c55e3fSMatthew G. Knepley     DM                 dm;
9418c55e3fSMatthew G. Knepley     PetscDS            prob;
9518c55e3fSMatthew G. Knepley     PetscSection       s;
967e01ee02SMatthew G. Knepley     Vec                Xstar = NULL, G = NULL;
9718c55e3fSMatthew G. Knepley     const PetscScalar *ax;
9818c55e3fSMatthew G. Knepley     PetscScalar       *axstar;
9918c55e3fSMatthew G. Knepley     PetscInt           Nf, f, pStart, pEnd, p;
10018c55e3fSMatthew G. Knepley 
10118c55e3fSMatthew G. Knepley     ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
10218c55e3fSMatthew G. Knepley     ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
10318c55e3fSMatthew G. Knepley     ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
10418c55e3fSMatthew G. Knepley     ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
10518c55e3fSMatthew G. Knepley     ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
10618c55e3fSMatthew G. Knepley     ierr = TSMimexGetXstarAndG(ts, dm, &Xstar, &G);CHKERRQ(ierr);
10718c55e3fSMatthew G. Knepley     ierr = VecCopy(X0, Xstar);CHKERRQ(ierr);
10818c55e3fSMatthew G. Knepley     ierr = VecGetArrayRead(x, &ax);CHKERRQ(ierr);
10918c55e3fSMatthew G. Knepley     ierr = VecGetArray(Xstar, &axstar);CHKERRQ(ierr);
11018c55e3fSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
11118c55e3fSMatthew G. Knepley       PetscBool implicit;
11218c55e3fSMatthew G. Knepley 
11318c55e3fSMatthew G. Knepley       ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr);
11418c55e3fSMatthew G. Knepley       if (!implicit) continue;
11518c55e3fSMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
11618c55e3fSMatthew G. Knepley         PetscScalar *a, *axs;
11718c55e3fSMatthew G. Knepley         PetscInt     fdof, fcdof, d;
11818c55e3fSMatthew G. Knepley 
11918c55e3fSMatthew G. Knepley         ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
12018c55e3fSMatthew G. Knepley         ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr);
12118c55e3fSMatthew G. Knepley         ierr = DMPlexPointGlobalFieldRead(dm, p, f, ax, &a);CHKERRQ(ierr);
12218c55e3fSMatthew G. Knepley         ierr = DMPlexPointGlobalFieldRef(dm, p, f, axstar, &axs);CHKERRQ(ierr);
12318c55e3fSMatthew G. Knepley         for (d = 0; d < fdof-fcdof; ++d) axs[d] = a[d];
12418c55e3fSMatthew G. Knepley       }
12518c55e3fSMatthew G. Knepley     }
12618c55e3fSMatthew G. Knepley     ierr = VecRestoreArrayRead(x, &ax);CHKERRQ(ierr);
12718c55e3fSMatthew G. Knepley     ierr = VecRestoreArray(Xstar, &axstar);CHKERRQ(ierr);
12818c55e3fSMatthew G. Knepley     ierr = TSComputeRHSFunction(ts, ts->ptime, Xstar, G);CHKERRQ(ierr);
12918c55e3fSMatthew G. Knepley     ierr = VecAXPY(y, -1.0, G);CHKERRQ(ierr);
13018c55e3fSMatthew G. Knepley     ierr = TSMimexRestoreXstarAndG(ts, dm, &Xstar, &G);CHKERRQ(ierr);
13118c55e3fSMatthew G. Knepley   }
132abadfa0fSMatthew G. Knepley   ts->dm = dmsave;
133abadfa0fSMatthew G. Knepley   ierr   = TSMimexRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
134abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
135abadfa0fSMatthew G. Knepley }
136abadfa0fSMatthew G. Knepley 
137abadfa0fSMatthew G. Knepley #undef __FUNCT__
138abadfa0fSMatthew G. Knepley #define __FUNCT__ "SNESTSFormJacobian_Mimex"
139abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_Mimex(SNES snes, Vec x, Mat A, Mat B, TS ts)
140abadfa0fSMatthew G. Knepley {
141abadfa0fSMatthew G. Knepley   TS_Mimex      *mimex = (TS_Mimex *) ts->data;
142abadfa0fSMatthew G. Knepley   DM             dm, dmsave;
143abadfa0fSMatthew G. Knepley   Vec            Xdot;
144abadfa0fSMatthew G. Knepley   PetscReal      shift = 1./ts->time_step;
145abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
146abadfa0fSMatthew G. Knepley 
147abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
148abadfa0fSMatthew G. Knepley   /* th->Xdot has already been computed in SNESTSFormFunction_Mimex (SNES guarantees this) */
149abadfa0fSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
150abadfa0fSMatthew G. Knepley   ierr = TSMimexGetX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr);
151abadfa0fSMatthew G. Knepley 
152abadfa0fSMatthew G. Knepley   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
153abadfa0fSMatthew G. Knepley   dmsave = ts->dm;
154abadfa0fSMatthew G. Knepley   ts->dm = dm;
155abadfa0fSMatthew G. Knepley   ierr   = TSComputeIJacobian(ts, mimex->stage_time, x, Xdot, shift, A, B, PETSC_TRUE);CHKERRQ(ierr);
156abadfa0fSMatthew G. Knepley   ts->dm = dmsave;
157abadfa0fSMatthew G. Knepley   ierr   = TSMimexRestoreX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr);
158abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
159abadfa0fSMatthew G. Knepley }
160abadfa0fSMatthew G. Knepley 
161abadfa0fSMatthew G. Knepley #undef __FUNCT__
16218c55e3fSMatthew G. Knepley #define __FUNCT__ "TSStep_Mimex_Split"
16318c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Split(TS ts)
164abadfa0fSMatthew G. Knepley {
165abadfa0fSMatthew G. Knepley   TS_Mimex          *mimex = (TS_Mimex *) ts->data;
166abadfa0fSMatthew G. Knepley   DM                 dm;
167abadfa0fSMatthew G. Knepley   PetscDS            prob;
168abadfa0fSMatthew G. Knepley   PetscSection       s;
169abadfa0fSMatthew G. Knepley   Vec                sol = ts->vec_sol, update = mimex->update;
170abadfa0fSMatthew G. Knepley   const PetscScalar *aupdate;
171abadfa0fSMatthew G. Knepley   PetscScalar       *asol, dt = ts->time_step;
172abadfa0fSMatthew G. Knepley   PetscInt           Nf, f, pStart, pEnd, p;
173abadfa0fSMatthew G. Knepley   PetscErrorCode     ierr;
174abadfa0fSMatthew G. Knepley 
175abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
176abadfa0fSMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
177abadfa0fSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
178abadfa0fSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
179abadfa0fSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
180abadfa0fSMatthew G. Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);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) {
18918c55e3fSMatthew G. Knepley     PetscBool implicit;
190abadfa0fSMatthew G. Knepley 
19118c55e3fSMatthew G. Knepley     ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr);
19218c55e3fSMatthew 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) {
21118c55e3fSMatthew G. Knepley     PetscBool implicit;
212abadfa0fSMatthew G. Knepley 
21318c55e3fSMatthew G. Knepley     ierr = PetscDSGetImplicit(prob, f, &implicit);CHKERRQ(ierr);
21418c55e3fSMatthew 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 
23418c55e3fSMatthew G. Knepley 
23518c55e3fSMatthew G. Knepley #undef __FUNCT__
23618c55e3fSMatthew G. Knepley #define __FUNCT__ "TSStep_Mimex_Implicit"
23718c55e3fSMatthew G. Knepley /* Evalute F at U and G at U0 for explicit fields and U for implicit fields */
23818c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex_Implicit(TS ts)
23918c55e3fSMatthew G. Knepley {
24018c55e3fSMatthew G. Knepley   TS_Mimex      *mimex  = (TS_Mimex *) ts->data;
24118c55e3fSMatthew G. Knepley   Vec            sol    = ts->vec_sol;
24218c55e3fSMatthew G. Knepley   Vec            update = mimex->update;
24318c55e3fSMatthew G. Knepley   PetscErrorCode ierr;
24418c55e3fSMatthew G. Knepley 
24518c55e3fSMatthew G. Knepley   PetscFunctionBegin;
24618c55e3fSMatthew G. Knepley   ierr = TSPreStage(ts, ts->ptime);CHKERRQ(ierr);
24718c55e3fSMatthew G. Knepley   /* Compute implicit update */
24818c55e3fSMatthew G. Knepley   mimex->stage_time = ts->ptime + ts->time_step;
24918c55e3fSMatthew G. Knepley   ts->ptime += ts->time_step;
25018c55e3fSMatthew G. Knepley   ierr = VecCopy(sol, update);CHKERRQ(ierr);
25118c55e3fSMatthew G. Knepley   ierr = SNESSolve(ts->snes, NULL, update);CHKERRQ(ierr);
25218c55e3fSMatthew G. Knepley   ierr = VecCopy(update, sol);CHKERRQ(ierr);
25318c55e3fSMatthew G. Knepley   ierr = TSPostStage(ts, ts->ptime, 0, &sol);CHKERRQ(ierr);
25418c55e3fSMatthew G. Knepley   ts->steps++;
25518c55e3fSMatthew G. Knepley   PetscFunctionReturn(0);
25618c55e3fSMatthew G. Knepley }
25718c55e3fSMatthew G. Knepley 
25818c55e3fSMatthew G. Knepley #undef __FUNCT__
25918c55e3fSMatthew G. Knepley #define __FUNCT__ "TSStep_Mimex"
26018c55e3fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex(TS ts)
26118c55e3fSMatthew G. Knepley {
26218c55e3fSMatthew G. Knepley   TS_Mimex       *mimex = (TS_Mimex*)ts->data;
26318c55e3fSMatthew G. Knepley   PetscErrorCode  ierr;
26418c55e3fSMatthew G. Knepley 
26518c55e3fSMatthew G. Knepley   PetscFunctionBegin;
26618c55e3fSMatthew G. Knepley   switch(mimex->version) {
26718c55e3fSMatthew G. Knepley   case 0:
26818c55e3fSMatthew G. Knepley     ierr = TSStep_Mimex_Split(ts);CHKERRQ(ierr); break;
26918c55e3fSMatthew G. Knepley   case 1:
27018c55e3fSMatthew G. Knepley     ierr = TSStep_Mimex_Implicit(ts);CHKERRQ(ierr); break;
27118c55e3fSMatthew G. Knepley   default:
27218c55e3fSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_OUTOFRANGE, "Unknown MIMEX version %d", mimex->version);
27318c55e3fSMatthew G. Knepley   }
27418c55e3fSMatthew G. Knepley   PetscFunctionReturn(0);
27518c55e3fSMatthew G. Knepley }
27618c55e3fSMatthew G. Knepley 
277abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/
278abadfa0fSMatthew G. Knepley 
279abadfa0fSMatthew G. Knepley #undef __FUNCT__
280abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSSetUp_Mimex"
281abadfa0fSMatthew G. Knepley static PetscErrorCode TSSetUp_Mimex(TS ts)
282abadfa0fSMatthew G. Knepley {
283abadfa0fSMatthew G. Knepley   TS_Mimex       *mimex = (TS_Mimex*)ts->data;
284abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
285abadfa0fSMatthew G. Knepley 
286abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
287abadfa0fSMatthew G. Knepley   ierr = VecDuplicate(ts->vec_sol, &mimex->update);CHKERRQ(ierr);
288abadfa0fSMatthew G. Knepley   ierr = VecDuplicate(ts->vec_sol, &mimex->Xdot);CHKERRQ(ierr);
289abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
290abadfa0fSMatthew G. Knepley }
291abadfa0fSMatthew G. Knepley 
292abadfa0fSMatthew G. Knepley #undef __FUNCT__
293abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSReset_Mimex"
294abadfa0fSMatthew G. Knepley static PetscErrorCode TSReset_Mimex(TS ts)
295abadfa0fSMatthew G. Knepley {
296abadfa0fSMatthew G. Knepley   TS_Mimex       *mimex = (TS_Mimex*)ts->data;
297abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
298abadfa0fSMatthew G. Knepley 
299abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
300abadfa0fSMatthew G. Knepley   ierr = VecDestroy(&mimex->update);CHKERRQ(ierr);
301abadfa0fSMatthew G. Knepley   ierr = VecDestroy(&mimex->Xdot);CHKERRQ(ierr);
302abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
303abadfa0fSMatthew G. Knepley }
304abadfa0fSMatthew G. Knepley 
305abadfa0fSMatthew G. Knepley #undef __FUNCT__
306abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSDestroy_Mimex"
307abadfa0fSMatthew G. Knepley static PetscErrorCode TSDestroy_Mimex(TS ts)
308abadfa0fSMatthew G. Knepley {
309abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
310abadfa0fSMatthew G. Knepley 
311abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
312abadfa0fSMatthew G. Knepley   ierr = TSReset_Mimex(ts);CHKERRQ(ierr);
313abadfa0fSMatthew G. Knepley   ierr = PetscFree(ts->data);CHKERRQ(ierr);
314abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
315abadfa0fSMatthew G. Knepley }
316abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/
317abadfa0fSMatthew G. Knepley 
318abadfa0fSMatthew G. Knepley #undef __FUNCT__
319abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSSetFromOptions_Mimex"
3204416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Mimex(PetscOptionItems *PetscOptionsObject, TS ts)
321abadfa0fSMatthew G. Knepley {
32218c55e3fSMatthew G. Knepley   TS_Mimex      *mimex = (TS_Mimex *) ts->data;
32318c55e3fSMatthew G. Knepley   PetscErrorCode ierr;
32418c55e3fSMatthew G. Knepley 
325abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
32618c55e3fSMatthew G. Knepley   ierr = PetscOptionsHead(PetscOptionsObject, "MIMEX ODE solver options");CHKERRQ(ierr);
32718c55e3fSMatthew G. Knepley   {
32818c55e3fSMatthew G. Knepley     ierr = PetscOptionsInt("-ts_mimex_version", "Algorithm version", "TSMimexSetVersion", mimex->version, &mimex->version, NULL);CHKERRQ(ierr);
32918c55e3fSMatthew G. Knepley   }
33018c55e3fSMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
331abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
332abadfa0fSMatthew G. Knepley }
333abadfa0fSMatthew G. Knepley 
334abadfa0fSMatthew G. Knepley #undef __FUNCT__
335abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSView_Mimex"
336abadfa0fSMatthew G. Knepley static PetscErrorCode TSView_Mimex(TS ts,PetscViewer viewer)
337abadfa0fSMatthew G. Knepley {
3385bf398f3SMatthew G. Knepley   TS_Mimex      *mimex = (TS_Mimex *) ts->data;
3395bf398f3SMatthew G. Knepley   PetscBool      iascii;
3405bf398f3SMatthew G. Knepley   PetscErrorCode ierr;
3415bf398f3SMatthew G. Knepley 
342abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
3435bf398f3SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
3445bf398f3SMatthew G. Knepley   if (iascii) {
3455bf398f3SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "  Version = %D\n", mimex->version);CHKERRQ(ierr);
3465bf398f3SMatthew G. Knepley   }
3475bf398f3SMatthew G. Knepley   if (ts->snes) {ierr = SNESView(ts->snes, viewer);CHKERRQ(ierr);}
348abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
349abadfa0fSMatthew G. Knepley }
350abadfa0fSMatthew G. Knepley 
351abadfa0fSMatthew G. Knepley #undef __FUNCT__
352abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSInterpolate_Mimex"
353abadfa0fSMatthew G. Knepley static PetscErrorCode TSInterpolate_Mimex(TS ts,PetscReal t,Vec X)
354abadfa0fSMatthew G. Knepley {
355abadfa0fSMatthew G. Knepley   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
356abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
357abadfa0fSMatthew G. Knepley 
358abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
359abadfa0fSMatthew G. Knepley   ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr);
360abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
361abadfa0fSMatthew G. Knepley }
362abadfa0fSMatthew G. Knepley 
363abadfa0fSMatthew G. Knepley #undef __FUNCT__
364abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSComputeLinearStability_Mimex"
365*560360afSLisandro Dalcin static PetscErrorCode TSComputeLinearStability_Mimex(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
366abadfa0fSMatthew G. Knepley {
367abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
368abadfa0fSMatthew G. Knepley   *yr = 1.0 + xr;
369abadfa0fSMatthew G. Knepley   *yi = xi;
370abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
371abadfa0fSMatthew G. Knepley }
372abadfa0fSMatthew G. Knepley /* ------------------------------------------------------------ */
373abadfa0fSMatthew G. Knepley 
374abadfa0fSMatthew G. Knepley /*MC
375abadfa0fSMatthew G. Knepley       TSMIMEX - ODE solver using the explicit forward Mimex method
376abadfa0fSMatthew G. Knepley 
377abadfa0fSMatthew G. Knepley   Level: beginner
378abadfa0fSMatthew G. Knepley 
379abadfa0fSMatthew G. Knepley .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
380abadfa0fSMatthew G. Knepley 
381abadfa0fSMatthew G. Knepley M*/
382abadfa0fSMatthew G. Knepley #undef __FUNCT__
383abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSCreate_Mimex"
384abadfa0fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_Mimex(TS ts)
385abadfa0fSMatthew G. Knepley {
386abadfa0fSMatthew G. Knepley   TS_Mimex       *mimex;
387abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
388abadfa0fSMatthew G. Knepley 
389abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
390abadfa0fSMatthew G. Knepley   ts->ops->setup           = TSSetUp_Mimex;
391abadfa0fSMatthew G. Knepley   ts->ops->step            = TSStep_Mimex;
392abadfa0fSMatthew G. Knepley   ts->ops->reset           = TSReset_Mimex;
393abadfa0fSMatthew G. Knepley   ts->ops->destroy         = TSDestroy_Mimex;
394abadfa0fSMatthew G. Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Mimex;
395abadfa0fSMatthew G. Knepley   ts->ops->view            = TSView_Mimex;
396abadfa0fSMatthew G. Knepley   ts->ops->interpolate     = TSInterpolate_Mimex;
397abadfa0fSMatthew G. Knepley   ts->ops->linearstability = TSComputeLinearStability_Mimex;
398abadfa0fSMatthew G. Knepley   ts->ops->snesfunction    = SNESTSFormFunction_Mimex;
399abadfa0fSMatthew G. Knepley   ts->ops->snesjacobian    = SNESTSFormJacobian_Mimex;
400abadfa0fSMatthew G. Knepley 
401abadfa0fSMatthew G. Knepley   ierr = PetscNewLog(ts,&mimex);CHKERRQ(ierr);
402abadfa0fSMatthew G. Knepley   ts->data = (void*)mimex;
403fdc0a3ceSMatthew G. Knepley 
404fdc0a3ceSMatthew G. Knepley   mimex->version = 1;
405abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
406abadfa0fSMatthew G. Knepley }
407