xref: /petsc/src/ts/impls/mimex/mimex.c (revision abadfa0ff52f084808973d56e72ca5b14b37a6b3)
1*abadfa0fSMatthew G. Knepley /*
2*abadfa0fSMatthew G. Knepley        Code for Timestepping with my makeshift IMEX.
3*abadfa0fSMatthew G. Knepley */
4*abadfa0fSMatthew G. Knepley #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
5*abadfa0fSMatthew G. Knepley #include <petscds.h>
6*abadfa0fSMatthew G. Knepley #include <petscdmplex.h>
7*abadfa0fSMatthew G. Knepley 
8*abadfa0fSMatthew G. Knepley typedef struct {
9*abadfa0fSMatthew G. Knepley   Vec       Xdot, update;
10*abadfa0fSMatthew G. Knepley   PetscReal stage_time;
11*abadfa0fSMatthew G. Knepley } TS_Mimex;
12*abadfa0fSMatthew G. Knepley 
13*abadfa0fSMatthew G. Knepley #undef __FUNCT__
14*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSMimexGetX0AndXdot"
15*abadfa0fSMatthew G. Knepley static PetscErrorCode TSMimexGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
16*abadfa0fSMatthew G. Knepley {
17*abadfa0fSMatthew G. Knepley   TS_Mimex      *mimex = (TS_Mimex *) ts->data;
18*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
19*abadfa0fSMatthew G. Knepley 
20*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
21*abadfa0fSMatthew G. Knepley   if (X0) {
22*abadfa0fSMatthew G. Knepley     if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSMimex_X0", X0);CHKERRQ(ierr);}
23*abadfa0fSMatthew G. Knepley     else                    {*X0  = ts->vec_sol;}
24*abadfa0fSMatthew G. Knepley   }
25*abadfa0fSMatthew G. Knepley   if (Xdot) {
26*abadfa0fSMatthew G. Knepley     if (dm && dm != ts->dm) {ierr  = DMGetNamedGlobalVector(dm, "TSMimex_Xdot", Xdot);CHKERRQ(ierr);}
27*abadfa0fSMatthew G. Knepley     else                    {*Xdot = mimex->Xdot;}
28*abadfa0fSMatthew G. Knepley   }
29*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
30*abadfa0fSMatthew G. Knepley }
31*abadfa0fSMatthew G. Knepley 
32*abadfa0fSMatthew G. Knepley 
33*abadfa0fSMatthew G. Knepley #undef __FUNCT__
34*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSMimexRestoreX0AndXdot"
35*abadfa0fSMatthew G. Knepley static PetscErrorCode TSMimexRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
36*abadfa0fSMatthew G. Knepley {
37*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
38*abadfa0fSMatthew G. Knepley 
39*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
40*abadfa0fSMatthew G. Knepley   if (X0)   if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_X0", X0);CHKERRQ(ierr);}
41*abadfa0fSMatthew G. Knepley   if (Xdot) if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSMimex_Xdot", Xdot);CHKERRQ(ierr);}
42*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
43*abadfa0fSMatthew G. Knepley }
44*abadfa0fSMatthew G. Knepley 
45*abadfa0fSMatthew G. Knepley /*
46*abadfa0fSMatthew G. Knepley   This defines the nonlinear equation that is to be solved with SNES
47*abadfa0fSMatthew G. Knepley   G(U) = F[t0+dt, U, (U-U0)*shift] = 0
48*abadfa0fSMatthew G. Knepley */
49*abadfa0fSMatthew G. Knepley #undef __FUNCT__
50*abadfa0fSMatthew G. Knepley #define __FUNCT__ "SNESTSFormFunction_Mimex"
51*abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_Mimex(SNES snes, Vec x, Vec y, TS ts)
52*abadfa0fSMatthew G. Knepley {
53*abadfa0fSMatthew G. Knepley   TS_Mimex      *mimex = (TS_Mimex *) ts->data;
54*abadfa0fSMatthew G. Knepley   DM             dm, dmsave;
55*abadfa0fSMatthew G. Knepley   Vec            X0, Xdot;
56*abadfa0fSMatthew G. Knepley   PetscReal      shift = 1./ts->time_step;
57*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
58*abadfa0fSMatthew G. Knepley 
59*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
60*abadfa0fSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
61*abadfa0fSMatthew G. Knepley   ierr = TSMimexGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
62*abadfa0fSMatthew G. Knepley   ierr = VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x);CHKERRQ(ierr);
63*abadfa0fSMatthew G. Knepley 
64*abadfa0fSMatthew G. Knepley   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
65*abadfa0fSMatthew G. Knepley   dmsave = ts->dm;
66*abadfa0fSMatthew G. Knepley   ts->dm = dm;
67*abadfa0fSMatthew G. Knepley   ierr   = TSComputeIFunction(ts, mimex->stage_time, x, Xdot, y, PETSC_TRUE);CHKERRQ(ierr);
68*abadfa0fSMatthew G. Knepley   ts->dm = dmsave;
69*abadfa0fSMatthew G. Knepley   ierr   = TSMimexRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
70*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
71*abadfa0fSMatthew G. Knepley }
72*abadfa0fSMatthew G. Knepley 
73*abadfa0fSMatthew G. Knepley #undef __FUNCT__
74*abadfa0fSMatthew G. Knepley #define __FUNCT__ "SNESTSFormJacobian_Mimex"
75*abadfa0fSMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_Mimex(SNES snes, Vec x, Mat A, Mat B, TS ts)
76*abadfa0fSMatthew G. Knepley {
77*abadfa0fSMatthew G. Knepley   TS_Mimex      *mimex = (TS_Mimex *) ts->data;
78*abadfa0fSMatthew G. Knepley   DM             dm, dmsave;
79*abadfa0fSMatthew G. Knepley   Vec            Xdot;
80*abadfa0fSMatthew G. Knepley   PetscReal      shift = 1./ts->time_step;
81*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
82*abadfa0fSMatthew G. Knepley 
83*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
84*abadfa0fSMatthew G. Knepley   /* th->Xdot has already been computed in SNESTSFormFunction_Mimex (SNES guarantees this) */
85*abadfa0fSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
86*abadfa0fSMatthew G. Knepley   ierr = TSMimexGetX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr);
87*abadfa0fSMatthew G. Knepley 
88*abadfa0fSMatthew G. Knepley   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
89*abadfa0fSMatthew G. Knepley   dmsave = ts->dm;
90*abadfa0fSMatthew G. Knepley   ts->dm = dm;
91*abadfa0fSMatthew G. Knepley   ierr   = TSComputeIJacobian(ts, mimex->stage_time, x, Xdot, shift, A, B, PETSC_TRUE);CHKERRQ(ierr);
92*abadfa0fSMatthew G. Knepley   ts->dm = dmsave;
93*abadfa0fSMatthew G. Knepley   ierr   = TSMimexRestoreX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr);
94*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
95*abadfa0fSMatthew G. Knepley }
96*abadfa0fSMatthew G. Knepley 
97*abadfa0fSMatthew G. Knepley #undef __FUNCT__
98*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSStep_Mimex"
99*abadfa0fSMatthew G. Knepley static PetscErrorCode TSStep_Mimex(TS ts)
100*abadfa0fSMatthew G. Knepley {
101*abadfa0fSMatthew G. Knepley   TS_Mimex          *mimex = (TS_Mimex *) ts->data;
102*abadfa0fSMatthew G. Knepley   DM                 dm;
103*abadfa0fSMatthew G. Knepley   PetscDS            prob;
104*abadfa0fSMatthew G. Knepley   PetscSection       s;
105*abadfa0fSMatthew G. Knepley   Vec                sol = ts->vec_sol, update = mimex->update;
106*abadfa0fSMatthew G. Knepley   const PetscScalar *aupdate;
107*abadfa0fSMatthew G. Knepley   PetscScalar       *asol, dt = ts->time_step;
108*abadfa0fSMatthew G. Knepley   PetscInt           Nf, f, pStart, pEnd, p;
109*abadfa0fSMatthew G. Knepley   PetscErrorCode     ierr;
110*abadfa0fSMatthew G. Knepley 
111*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
112*abadfa0fSMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
113*abadfa0fSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
114*abadfa0fSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
115*abadfa0fSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
116*abadfa0fSMatthew G. Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
117*abadfa0fSMatthew G. Knepley   ierr = TSPreStep(ts);CHKERRQ(ierr);
118*abadfa0fSMatthew G. Knepley   ierr = TSPreStage(ts, ts->ptime);CHKERRQ(ierr);
119*abadfa0fSMatthew G. Knepley   /* Compute implicit update */
120*abadfa0fSMatthew G. Knepley   mimex->stage_time = ts->ptime + ts->time_step;
121*abadfa0fSMatthew G. Knepley   ierr = VecCopy(sol, update);CHKERRQ(ierr);
122*abadfa0fSMatthew G. Knepley   ierr = SNESSolve(ts->snes, NULL, update);CHKERRQ(ierr);
123*abadfa0fSMatthew G. Knepley   ierr = VecGetArrayRead(update, &aupdate);CHKERRQ(ierr);
124*abadfa0fSMatthew G. Knepley   ierr = VecGetArray(sol, &asol);CHKERRQ(ierr);
125*abadfa0fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
126*abadfa0fSMatthew G. Knepley     PetscObject  obj;
127*abadfa0fSMatthew G. Knepley     PetscClassId id;
128*abadfa0fSMatthew G. Knepley     PetscInt     pdim, d;
129*abadfa0fSMatthew G. Knepley 
130*abadfa0fSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
131*abadfa0fSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
132*abadfa0fSMatthew G. Knepley     if (id == PETSCFV_CLASSID) continue;
133*abadfa0fSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
134*abadfa0fSMatthew G. Knepley       PetscScalar *au, *as;
135*abadfa0fSMatthew G. Knepley       PetscInt     fdof, fcdof, d;
136*abadfa0fSMatthew G. Knepley 
137*abadfa0fSMatthew G. Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
138*abadfa0fSMatthew G. Knepley       ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr);
139*abadfa0fSMatthew G. Knepley       ierr = DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au);CHKERRQ(ierr);
140*abadfa0fSMatthew G. Knepley       ierr = DMPlexPointGlobalFieldRef(dm, p, f, asol, &as);CHKERRQ(ierr);
141*abadfa0fSMatthew G. Knepley       for (d = 0; d < fdof-fcdof; ++d) as[d] = au[d];
142*abadfa0fSMatthew G. Knepley     }
143*abadfa0fSMatthew G. Knepley   }
144*abadfa0fSMatthew G. Knepley   ierr = VecRestoreArrayRead(update, &aupdate);CHKERRQ(ierr);
145*abadfa0fSMatthew G. Knepley   ierr = VecRestoreArray(sol, &asol);CHKERRQ(ierr);
146*abadfa0fSMatthew G. Knepley   /* Compute explicit update */
147*abadfa0fSMatthew G. Knepley   ierr = TSComputeRHSFunction(ts, ts->ptime, sol, update);CHKERRQ(ierr);
148*abadfa0fSMatthew G. Knepley   ierr = VecGetArrayRead(update, &aupdate);CHKERRQ(ierr);
149*abadfa0fSMatthew G. Knepley   ierr = VecGetArray(sol, &asol);CHKERRQ(ierr);
150*abadfa0fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
151*abadfa0fSMatthew G. Knepley     PetscObject  obj;
152*abadfa0fSMatthew G. Knepley     PetscClassId id;
153*abadfa0fSMatthew G. Knepley     PetscInt     pdim, d;
154*abadfa0fSMatthew G. Knepley 
155*abadfa0fSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
156*abadfa0fSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
157*abadfa0fSMatthew G. Knepley     if (id != PETSCFV_CLASSID) continue;
158*abadfa0fSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
159*abadfa0fSMatthew G. Knepley       PetscScalar *au, *as;
160*abadfa0fSMatthew G. Knepley       PetscInt     fdof, fcdof, d;
161*abadfa0fSMatthew G. Knepley 
162*abadfa0fSMatthew G. Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
163*abadfa0fSMatthew G. Knepley       ierr = PetscSectionGetFieldConstraintDof(s, p, f, &fcdof);CHKERRQ(ierr);
164*abadfa0fSMatthew G. Knepley       ierr = DMPlexPointGlobalFieldRead(dm, p, f, aupdate, &au);CHKERRQ(ierr);
165*abadfa0fSMatthew G. Knepley       ierr = DMPlexPointGlobalFieldRef(dm, p, f, asol, &as);CHKERRQ(ierr);
166*abadfa0fSMatthew G. Knepley       for (d = 0; d < fdof-fcdof; ++d) as[d] += dt*au[d];
167*abadfa0fSMatthew G. Knepley     }
168*abadfa0fSMatthew G. Knepley   }
169*abadfa0fSMatthew G. Knepley   ierr = VecRestoreArrayRead(update, &aupdate);CHKERRQ(ierr);
170*abadfa0fSMatthew G. Knepley   ierr = VecRestoreArray(sol, &asol);CHKERRQ(ierr);
171*abadfa0fSMatthew G. Knepley   ierr = TSPostStage(ts, ts->ptime, 0, &sol);CHKERRQ(ierr);
172*abadfa0fSMatthew G. Knepley   ts->ptime += ts->time_step;
173*abadfa0fSMatthew G. Knepley   ts->steps++;
174*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
175*abadfa0fSMatthew G. Knepley }
176*abadfa0fSMatthew G. Knepley 
177*abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/
178*abadfa0fSMatthew G. Knepley 
179*abadfa0fSMatthew G. Knepley #undef __FUNCT__
180*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSSetUp_Mimex"
181*abadfa0fSMatthew G. Knepley static PetscErrorCode TSSetUp_Mimex(TS ts)
182*abadfa0fSMatthew G. Knepley {
183*abadfa0fSMatthew G. Knepley   TS_Mimex       *mimex = (TS_Mimex*)ts->data;
184*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
185*abadfa0fSMatthew G. Knepley 
186*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
187*abadfa0fSMatthew G. Knepley   ierr = VecDuplicate(ts->vec_sol, &mimex->update);CHKERRQ(ierr);
188*abadfa0fSMatthew G. Knepley   ierr = VecDuplicate(ts->vec_sol, &mimex->Xdot);CHKERRQ(ierr);
189*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
190*abadfa0fSMatthew G. Knepley }
191*abadfa0fSMatthew G. Knepley 
192*abadfa0fSMatthew G. Knepley #undef __FUNCT__
193*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSReset_Mimex"
194*abadfa0fSMatthew G. Knepley static PetscErrorCode TSReset_Mimex(TS ts)
195*abadfa0fSMatthew G. Knepley {
196*abadfa0fSMatthew G. Knepley   TS_Mimex       *mimex = (TS_Mimex*)ts->data;
197*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
198*abadfa0fSMatthew G. Knepley 
199*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
200*abadfa0fSMatthew G. Knepley   ierr = VecDestroy(&mimex->update);CHKERRQ(ierr);
201*abadfa0fSMatthew G. Knepley   ierr = VecDestroy(&mimex->Xdot);CHKERRQ(ierr);
202*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
203*abadfa0fSMatthew G. Knepley }
204*abadfa0fSMatthew G. Knepley 
205*abadfa0fSMatthew G. Knepley #undef __FUNCT__
206*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSDestroy_Mimex"
207*abadfa0fSMatthew G. Knepley static PetscErrorCode TSDestroy_Mimex(TS ts)
208*abadfa0fSMatthew G. Knepley {
209*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
210*abadfa0fSMatthew G. Knepley 
211*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
212*abadfa0fSMatthew G. Knepley   ierr = TSReset_Mimex(ts);CHKERRQ(ierr);
213*abadfa0fSMatthew G. Knepley   ierr = PetscFree(ts->data);CHKERRQ(ierr);
214*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
215*abadfa0fSMatthew G. Knepley }
216*abadfa0fSMatthew G. Knepley /*------------------------------------------------------------*/
217*abadfa0fSMatthew G. Knepley 
218*abadfa0fSMatthew G. Knepley #undef __FUNCT__
219*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSSetFromOptions_Mimex"
220*abadfa0fSMatthew G. Knepley static PetscErrorCode TSSetFromOptions_Mimex(PetscOptions *PetscOptionsObject,TS ts)
221*abadfa0fSMatthew G. Knepley {
222*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
223*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
224*abadfa0fSMatthew G. Knepley }
225*abadfa0fSMatthew G. Knepley 
226*abadfa0fSMatthew G. Knepley #undef __FUNCT__
227*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSView_Mimex"
228*abadfa0fSMatthew G. Knepley static PetscErrorCode TSView_Mimex(TS ts,PetscViewer viewer)
229*abadfa0fSMatthew G. Knepley {
230*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
231*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
232*abadfa0fSMatthew G. Knepley }
233*abadfa0fSMatthew G. Knepley 
234*abadfa0fSMatthew G. Knepley #undef __FUNCT__
235*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSInterpolate_Mimex"
236*abadfa0fSMatthew G. Knepley static PetscErrorCode TSInterpolate_Mimex(TS ts,PetscReal t,Vec X)
237*abadfa0fSMatthew G. Knepley {
238*abadfa0fSMatthew G. Knepley   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
239*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
240*abadfa0fSMatthew G. Knepley 
241*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
242*abadfa0fSMatthew G. Knepley   ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr);
243*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
244*abadfa0fSMatthew G. Knepley }
245*abadfa0fSMatthew G. Knepley 
246*abadfa0fSMatthew G. Knepley #undef __FUNCT__
247*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSComputeLinearStability_Mimex"
248*abadfa0fSMatthew G. Knepley PetscErrorCode TSComputeLinearStability_Mimex(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
249*abadfa0fSMatthew G. Knepley {
250*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
251*abadfa0fSMatthew G. Knepley   *yr = 1.0 + xr;
252*abadfa0fSMatthew G. Knepley   *yi = xi;
253*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
254*abadfa0fSMatthew G. Knepley }
255*abadfa0fSMatthew G. Knepley /* ------------------------------------------------------------ */
256*abadfa0fSMatthew G. Knepley 
257*abadfa0fSMatthew G. Knepley /*MC
258*abadfa0fSMatthew G. Knepley       TSMIMEX - ODE solver using the explicit forward Mimex method
259*abadfa0fSMatthew G. Knepley 
260*abadfa0fSMatthew G. Knepley   Level: beginner
261*abadfa0fSMatthew G. Knepley 
262*abadfa0fSMatthew G. Knepley .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
263*abadfa0fSMatthew G. Knepley 
264*abadfa0fSMatthew G. Knepley M*/
265*abadfa0fSMatthew G. Knepley #undef __FUNCT__
266*abadfa0fSMatthew G. Knepley #define __FUNCT__ "TSCreate_Mimex"
267*abadfa0fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_Mimex(TS ts)
268*abadfa0fSMatthew G. Knepley {
269*abadfa0fSMatthew G. Knepley   TS_Mimex       *mimex;
270*abadfa0fSMatthew G. Knepley   PetscErrorCode ierr;
271*abadfa0fSMatthew G. Knepley 
272*abadfa0fSMatthew G. Knepley   PetscFunctionBegin;
273*abadfa0fSMatthew G. Knepley   ts->ops->setup           = TSSetUp_Mimex;
274*abadfa0fSMatthew G. Knepley   ts->ops->step            = TSStep_Mimex;
275*abadfa0fSMatthew G. Knepley   ts->ops->reset           = TSReset_Mimex;
276*abadfa0fSMatthew G. Knepley   ts->ops->destroy         = TSDestroy_Mimex;
277*abadfa0fSMatthew G. Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Mimex;
278*abadfa0fSMatthew G. Knepley   ts->ops->view            = TSView_Mimex;
279*abadfa0fSMatthew G. Knepley   ts->ops->interpolate     = TSInterpolate_Mimex;
280*abadfa0fSMatthew G. Knepley   ts->ops->linearstability = TSComputeLinearStability_Mimex;
281*abadfa0fSMatthew G. Knepley   ts->ops->snesfunction    = SNESTSFormFunction_Mimex;
282*abadfa0fSMatthew G. Knepley   ts->ops->snesjacobian    = SNESTSFormJacobian_Mimex;
283*abadfa0fSMatthew G. Knepley 
284*abadfa0fSMatthew G. Knepley   ierr = PetscNewLog(ts,&mimex);CHKERRQ(ierr);
285*abadfa0fSMatthew G. Knepley   ts->data = (void*)mimex;
286*abadfa0fSMatthew G. Knepley   PetscFunctionReturn(0);
287*abadfa0fSMatthew G. Knepley }
288