xref: /petsc/src/ts/impls/implicit/discgrad/tsdiscgrad.c (revision f1a722f8c411b335aaa5b04a3b4219fef8ff4c99)
140be0ff1SMatthew G. Knepley /*
240be0ff1SMatthew G. Knepley   Code for timestepping with discrete gradient integrators
340be0ff1SMatthew G. Knepley */
440be0ff1SMatthew G. Knepley #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
540be0ff1SMatthew G. Knepley #include <petscdm.h>
640be0ff1SMatthew G. Knepley 
740be0ff1SMatthew G. Knepley PetscBool DGCite = PETSC_FALSE;
840be0ff1SMatthew G. Knepley const char DGCitation[] = "@article{Gonzalez1996,\n"
940be0ff1SMatthew G. Knepley                           "  title   = {Time integration and discrete Hamiltonian systems},\n"
1040be0ff1SMatthew G. Knepley                           "  author  = {Oscar Gonzalez},\n"
1140be0ff1SMatthew G. Knepley                           "  journal = {Journal of Nonlinear Science},\n"
1240be0ff1SMatthew G. Knepley                           "  volume  = {6},\n"
1340be0ff1SMatthew G. Knepley                           "  pages   = {449--467},\n"
1440be0ff1SMatthew G. Knepley                           "  doi     = {10.1007/978-1-4612-1246-1_10},\n"
1540be0ff1SMatthew G. Knepley                           "  year    = {1996}\n}\n";
1640be0ff1SMatthew G. Knepley 
1740be0ff1SMatthew G. Knepley typedef struct {
1840be0ff1SMatthew G. Knepley   PetscReal    stage_time;
1940be0ff1SMatthew G. Knepley   Vec          X0, X, Xdot;
20db4653c2SDaniel Finn   void        *funcCtx;
21db4653c2SDaniel Finn   PetscBool    gonzalez;
2240be0ff1SMatthew G. Knepley   PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *);
2340be0ff1SMatthew G. Knepley   PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *);
2440be0ff1SMatthew G. Knepley   PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *);
2540be0ff1SMatthew G. Knepley } TS_DiscGrad;
2640be0ff1SMatthew G. Knepley 
2740be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGradGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
2840be0ff1SMatthew G. Knepley {
2940be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
3040be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
3140be0ff1SMatthew G. Knepley 
3240be0ff1SMatthew G. Knepley   PetscFunctionBegin;
3340be0ff1SMatthew G. Knepley   if (X0) {
3440be0ff1SMatthew G. Knepley     if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0);CHKERRQ(ierr);}
3540be0ff1SMatthew G. Knepley     else                    {*X0  = ts->vec_sol;}
3640be0ff1SMatthew G. Knepley   }
3740be0ff1SMatthew G. Knepley   if (Xdot) {
3840be0ff1SMatthew G. Knepley     if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot);CHKERRQ(ierr);}
3940be0ff1SMatthew G. Knepley     else                    {*Xdot = dg->Xdot;}
4040be0ff1SMatthew G. Knepley   }
4140be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
4240be0ff1SMatthew G. Knepley }
4340be0ff1SMatthew G. Knepley 
4440be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
4540be0ff1SMatthew G. Knepley {
4640be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
4740be0ff1SMatthew G. Knepley 
4840be0ff1SMatthew G. Knepley   PetscFunctionBegin;
4940be0ff1SMatthew G. Knepley   if (X0) {
5040be0ff1SMatthew G. Knepley     if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0);CHKERRQ(ierr);}
5140be0ff1SMatthew G. Knepley   }
5240be0ff1SMatthew G. Knepley   if (Xdot) {
5340be0ff1SMatthew G. Knepley     if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot);CHKERRQ(ierr);}
5440be0ff1SMatthew G. Knepley   }
5540be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
5640be0ff1SMatthew G. Knepley }
5740be0ff1SMatthew G. Knepley 
5840be0ff1SMatthew G. Knepley static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, void *ctx)
5940be0ff1SMatthew G. Knepley {
6040be0ff1SMatthew G. Knepley   PetscFunctionBegin;
6140be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
6240be0ff1SMatthew G. Knepley }
6340be0ff1SMatthew G. Knepley 
6440be0ff1SMatthew G. Knepley static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
6540be0ff1SMatthew G. Knepley {
6640be0ff1SMatthew G. Knepley   TS             ts = (TS) ctx;
6740be0ff1SMatthew G. Knepley   Vec            X0, Xdot, X0_c, Xdot_c;
6840be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
6940be0ff1SMatthew G. Knepley 
7040be0ff1SMatthew G. Knepley   PetscFunctionBegin;
7140be0ff1SMatthew G. Knepley   ierr = TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot);CHKERRQ(ierr);
7240be0ff1SMatthew G. Knepley   ierr = TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c);CHKERRQ(ierr);
7340be0ff1SMatthew G. Knepley   ierr = MatRestrict(restrct, X0, X0_c);CHKERRQ(ierr);
7440be0ff1SMatthew G. Knepley   ierr = MatRestrict(restrct, Xdot, Xdot_c);CHKERRQ(ierr);
7540be0ff1SMatthew G. Knepley   ierr = VecPointwiseMult(X0_c, rscale, X0_c);CHKERRQ(ierr);
7640be0ff1SMatthew G. Knepley   ierr = VecPointwiseMult(Xdot_c, rscale, Xdot_c);CHKERRQ(ierr);
7740be0ff1SMatthew G. Knepley   ierr = TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot);CHKERRQ(ierr);
7840be0ff1SMatthew G. Knepley   ierr = TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c);CHKERRQ(ierr);
7940be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
8040be0ff1SMatthew G. Knepley }
8140be0ff1SMatthew G. Knepley 
8240be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, void *ctx)
8340be0ff1SMatthew G. Knepley {
8440be0ff1SMatthew G. Knepley   PetscFunctionBegin;
8540be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
8640be0ff1SMatthew G. Knepley }
8740be0ff1SMatthew G. Knepley 
8840be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
8940be0ff1SMatthew G. Knepley {
9040be0ff1SMatthew G. Knepley   TS             ts = (TS) ctx;
9140be0ff1SMatthew G. Knepley   Vec            X0, Xdot, X0_sub, Xdot_sub;
9240be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
9340be0ff1SMatthew G. Knepley 
9440be0ff1SMatthew G. Knepley   PetscFunctionBegin;
9540be0ff1SMatthew G. Knepley   ierr = TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
9640be0ff1SMatthew G. Knepley   ierr = TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub);CHKERRQ(ierr);
9740be0ff1SMatthew G. Knepley 
9840be0ff1SMatthew G. Knepley   ierr = VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
9940be0ff1SMatthew G. Knepley   ierr = VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
10040be0ff1SMatthew G. Knepley 
10140be0ff1SMatthew G. Knepley   ierr = VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
10240be0ff1SMatthew G. Knepley   ierr = VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
10340be0ff1SMatthew G. Knepley 
10440be0ff1SMatthew G. Knepley   ierr = TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
10540be0ff1SMatthew G. Knepley   ierr = TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub);CHKERRQ(ierr);
10640be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
10740be0ff1SMatthew G. Knepley }
10840be0ff1SMatthew G. Knepley 
10940be0ff1SMatthew G. Knepley static PetscErrorCode TSSetUp_DiscGrad(TS ts)
11040be0ff1SMatthew G. Knepley {
11140be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
11240be0ff1SMatthew G. Knepley   DM             dm;
11340be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
11440be0ff1SMatthew G. Knepley 
11540be0ff1SMatthew G. Knepley   PetscFunctionBegin;
11640be0ff1SMatthew G. Knepley   if (!dg->X)    {ierr = VecDuplicate(ts->vec_sol, &dg->X);CHKERRQ(ierr);}
11740be0ff1SMatthew G. Knepley   if (!dg->X0)   {ierr = VecDuplicate(ts->vec_sol, &dg->X0);CHKERRQ(ierr);}
11840be0ff1SMatthew G. Knepley   if (!dg->Xdot) {ierr = VecDuplicate(ts->vec_sol, &dg->Xdot);CHKERRQ(ierr);}
11940be0ff1SMatthew G. Knepley 
12040be0ff1SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
12140be0ff1SMatthew G. Knepley   ierr = DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr);
12240be0ff1SMatthew G. Knepley   ierr = DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr);
12340be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
12440be0ff1SMatthew G. Knepley }
12540be0ff1SMatthew G. Knepley 
12640be0ff1SMatthew G. Knepley static PetscErrorCode TSSetFromOptions_DiscGrad(PetscOptionItems *PetscOptionsObject, TS ts)
12740be0ff1SMatthew G. Knepley {
128db4653c2SDaniel Finn   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
12940be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
13040be0ff1SMatthew G. Knepley 
13140be0ff1SMatthew G. Knepley   PetscFunctionBegin;
13240be0ff1SMatthew G. Knepley   ierr = PetscOptionsHead(PetscOptionsObject, "Discrete Gradients ODE solver options");CHKERRQ(ierr);
13340be0ff1SMatthew G. Knepley   {
134db4653c2SDaniel Finn     ierr = PetscOptionsBool("-ts_discgrad_gonzalez","Use Gonzalez term in discrete gradients formulation","TSDiscGradUseGonzalez",dg->gonzalez,&dg->gonzalez,NULL);CHKERRQ(ierr);
13540be0ff1SMatthew G. Knepley   }
13640be0ff1SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
13740be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
13840be0ff1SMatthew G. Knepley }
13940be0ff1SMatthew G. Knepley 
14040be0ff1SMatthew G. Knepley static PetscErrorCode TSView_DiscGrad(TS ts,PetscViewer viewer)
14140be0ff1SMatthew G. Knepley {
14240be0ff1SMatthew G. Knepley   PetscBool      iascii;
14340be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
14440be0ff1SMatthew G. Knepley 
14540be0ff1SMatthew G. Knepley   PetscFunctionBegin;
14640be0ff1SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
14740be0ff1SMatthew G. Knepley   if (iascii) {
14840be0ff1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer,"  Discrete Gradients\n");CHKERRQ(ierr);
14940be0ff1SMatthew G. Knepley   }
15040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
15140be0ff1SMatthew G. Knepley }
15240be0ff1SMatthew G. Knepley 
153db4653c2SDaniel Finn static PetscErrorCode TSDiscGradIsGonzalez_DiscGrad(TS ts,PetscBool *gonzalez)
154db4653c2SDaniel Finn {
155db4653c2SDaniel Finn   TS_DiscGrad *dg = (TS_DiscGrad*)ts->data;
156db4653c2SDaniel Finn 
157db4653c2SDaniel Finn   PetscFunctionBegin;
158db4653c2SDaniel Finn   *gonzalez = dg->gonzalez;
159db4653c2SDaniel Finn   PetscFunctionReturn(0);
160db4653c2SDaniel Finn }
161db4653c2SDaniel Finn 
162db4653c2SDaniel Finn static PetscErrorCode TSDiscGradUseGonzalez_DiscGrad(TS ts,PetscBool flg)
163db4653c2SDaniel Finn {
164db4653c2SDaniel Finn   TS_DiscGrad *dg = (TS_DiscGrad*)ts->data;
165db4653c2SDaniel Finn 
166db4653c2SDaniel Finn   PetscFunctionBegin;
167db4653c2SDaniel Finn   dg->gonzalez = flg;
168db4653c2SDaniel Finn   PetscFunctionReturn(0);
169db4653c2SDaniel Finn }
170db4653c2SDaniel Finn 
17140be0ff1SMatthew G. Knepley static PetscErrorCode TSReset_DiscGrad(TS ts)
17240be0ff1SMatthew G. Knepley {
17340be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
17440be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
17540be0ff1SMatthew G. Knepley 
17640be0ff1SMatthew G. Knepley   PetscFunctionBegin;
17740be0ff1SMatthew G. Knepley   ierr = VecDestroy(&dg->X);CHKERRQ(ierr);
17840be0ff1SMatthew G. Knepley   ierr = VecDestroy(&dg->X0);CHKERRQ(ierr);
17940be0ff1SMatthew G. Knepley   ierr = VecDestroy(&dg->Xdot);CHKERRQ(ierr);
18040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
18140be0ff1SMatthew G. Knepley }
18240be0ff1SMatthew G. Knepley 
18340be0ff1SMatthew G. Knepley static PetscErrorCode TSDestroy_DiscGrad(TS ts)
18440be0ff1SMatthew G. Knepley {
18540be0ff1SMatthew G. Knepley   DM             dm;
18640be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
18740be0ff1SMatthew G. Knepley 
18840be0ff1SMatthew G. Knepley   PetscFunctionBegin;
18940be0ff1SMatthew G. Knepley   ierr = TSReset_DiscGrad(ts);CHKERRQ(ierr);
19040be0ff1SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
19140be0ff1SMatthew G. Knepley   if (dm) {
19240be0ff1SMatthew G. Knepley     ierr = DMCoarsenHookRemove(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr);
19340be0ff1SMatthew G. Knepley     ierr = DMSubDomainHookRemove(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr);
19440be0ff1SMatthew G. Knepley   }
19540be0ff1SMatthew G. Knepley   ierr = PetscFree(ts->data);CHKERRQ(ierr);
19640be0ff1SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",NULL);CHKERRQ(ierr);
19740be0ff1SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",NULL);CHKERRQ(ierr);
19840be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
19940be0ff1SMatthew G. Knepley }
20040be0ff1SMatthew G. Knepley 
20140be0ff1SMatthew G. Knepley static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X)
20240be0ff1SMatthew G. Knepley {
20340be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad*)ts->data;
20440be0ff1SMatthew G. Knepley   PetscReal      dt = t - ts->ptime;
20540be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
20640be0ff1SMatthew G. Knepley 
20740be0ff1SMatthew G. Knepley   PetscFunctionBegin;
20840be0ff1SMatthew G. Knepley   ierr = VecCopy(ts->vec_sol, dg->X);CHKERRQ(ierr);
20940be0ff1SMatthew G. Knepley   ierr = VecWAXPY(X, dt, dg->Xdot, dg->X);CHKERRQ(ierr);
21040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
21140be0ff1SMatthew G. Knepley }
21240be0ff1SMatthew G. Knepley 
21340be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x)
21440be0ff1SMatthew G. Knepley {
21540be0ff1SMatthew G. Knepley   SNES           snes;
21640be0ff1SMatthew G. Knepley   PetscInt       nits, lits;
21740be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
21840be0ff1SMatthew G. Knepley 
21940be0ff1SMatthew G. Knepley   PetscFunctionBegin;
22040be0ff1SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
22140be0ff1SMatthew G. Knepley   ierr = SNESSolve(snes, b, x);CHKERRQ(ierr);
22240be0ff1SMatthew G. Knepley   ierr = SNESGetIterationNumber(snes, &nits);CHKERRQ(ierr);
22340be0ff1SMatthew G. Knepley   ierr = SNESGetLinearSolveIterations(snes, &lits);CHKERRQ(ierr);
22440be0ff1SMatthew G. Knepley   ts->snes_its += nits;
22540be0ff1SMatthew G. Knepley   ts->ksp_its  += lits;
22640be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
22740be0ff1SMatthew G. Knepley }
22840be0ff1SMatthew G. Knepley 
22940be0ff1SMatthew G. Knepley static PetscErrorCode TSStep_DiscGrad(TS ts)
23040be0ff1SMatthew G. Knepley {
23140be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
23240be0ff1SMatthew G. Knepley   TSAdapt        adapt;
23340be0ff1SMatthew G. Knepley   TSStepStatus   status          = TS_STEP_INCOMPLETE;
23440be0ff1SMatthew G. Knepley   PetscInt       rejections      = 0;
23540be0ff1SMatthew G. Knepley   PetscBool      stageok, accept = PETSC_TRUE;
23640be0ff1SMatthew G. Knepley   PetscReal      next_time_step  = ts->time_step;
23740be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
23840be0ff1SMatthew G. Knepley 
23940be0ff1SMatthew G. Knepley   PetscFunctionBegin;
24040be0ff1SMatthew G. Knepley   ierr = TSGetAdapt(ts, &adapt);CHKERRQ(ierr);
24140be0ff1SMatthew G. Knepley   if (!ts->steprollback) {ierr = VecCopy(ts->vec_sol, dg->X0);CHKERRQ(ierr);}
24240be0ff1SMatthew G. Knepley 
24340be0ff1SMatthew G. Knepley   while (!ts->reason && status != TS_STEP_COMPLETE) {
24440be0ff1SMatthew G. Knepley     PetscReal shift = 1/(0.5*ts->time_step);
24540be0ff1SMatthew G. Knepley 
24640be0ff1SMatthew G. Knepley     dg->stage_time = ts->ptime + 0.5*ts->time_step;
24740be0ff1SMatthew G. Knepley 
24840be0ff1SMatthew G. Knepley     ierr = VecCopy(dg->X0, dg->X);CHKERRQ(ierr);
24940be0ff1SMatthew G. Knepley     ierr = TSPreStage(ts, dg->stage_time);CHKERRQ(ierr);
25040be0ff1SMatthew G. Knepley     ierr = TSDiscGrad_SNESSolve(ts, NULL, dg->X);CHKERRQ(ierr);
25140be0ff1SMatthew G. Knepley     ierr = TSPostStage(ts, dg->stage_time, 0, &dg->X);CHKERRQ(ierr);
25240be0ff1SMatthew G. Knepley     ierr = TSAdaptCheckStage(adapt, ts, dg->stage_time, dg->X, &stageok);CHKERRQ(ierr);
25340be0ff1SMatthew G. Knepley     if (!stageok) goto reject_step;
25440be0ff1SMatthew G. Knepley 
25540be0ff1SMatthew G. Knepley     status = TS_STEP_PENDING;
25640be0ff1SMatthew G. Knepley     ierr = VecAXPBYPCZ(dg->Xdot, -shift, shift, 0, dg->X0, dg->X);CHKERRQ(ierr);
25740be0ff1SMatthew G. Knepley     ierr = VecAXPY(ts->vec_sol, ts->time_step, dg->Xdot);CHKERRQ(ierr);
25840be0ff1SMatthew G. Knepley     ierr = TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept);CHKERRQ(ierr);
25940be0ff1SMatthew G. Knepley     status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
26040be0ff1SMatthew G. Knepley     if (!accept) {
26140be0ff1SMatthew G. Knepley       ierr = VecCopy(dg->X0, ts->vec_sol);CHKERRQ(ierr);
26240be0ff1SMatthew G. Knepley       ts->time_step = next_time_step;
26340be0ff1SMatthew G. Knepley       goto reject_step;
26440be0ff1SMatthew G. Knepley     }
26540be0ff1SMatthew G. Knepley     ts->ptime    += ts->time_step;
26640be0ff1SMatthew G. Knepley     ts->time_step = next_time_step;
26740be0ff1SMatthew G. Knepley     break;
26840be0ff1SMatthew G. Knepley 
26940be0ff1SMatthew G. Knepley   reject_step:
27040be0ff1SMatthew G. Knepley     ts->reject++; accept = PETSC_FALSE;
27140be0ff1SMatthew G. Knepley     if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) {
27240be0ff1SMatthew G. Knepley       ts->reason = TS_DIVERGED_STEP_REJECTED;
2737d3de750SJacob Faibussowitsch       ierr = PetscInfo(ts, "Step=%D, step rejections %D greater than current TS allowed, stopping solve\n", ts->steps, rejections);CHKERRQ(ierr);
27440be0ff1SMatthew G. Knepley     }
27540be0ff1SMatthew G. Knepley   }
27640be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
27740be0ff1SMatthew G. Knepley }
27840be0ff1SMatthew G. Knepley 
27940be0ff1SMatthew G. Knepley static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y)
28040be0ff1SMatthew G. Knepley {
28140be0ff1SMatthew G. Knepley   TS_DiscGrad *dg = (TS_DiscGrad *) ts->data;
28240be0ff1SMatthew G. Knepley 
28340be0ff1SMatthew G. Knepley   PetscFunctionBegin;
28440be0ff1SMatthew G. Knepley   if (ns) *ns = 1;
28540be0ff1SMatthew G. Knepley   if (Y)  *Y  = &(dg->X);
28640be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
28740be0ff1SMatthew G. Knepley }
28840be0ff1SMatthew G. Knepley 
28940be0ff1SMatthew G. Knepley /*
29040be0ff1SMatthew G. Knepley   This defines the nonlinear equation that is to be solved with SNES
29140be0ff1SMatthew G. Knepley     G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0
29240be0ff1SMatthew G. Knepley */
293db4653c2SDaniel Finn 
294db4653c2SDaniel Finn /* x = (x+x')/2 */
295db4653c2SDaniel Finn /* NEED TO CALCULATE x_{n+1} from x and x_{n}*/
29640be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts)
29740be0ff1SMatthew G. Knepley {
298db4653c2SDaniel Finn 
29940be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg    = (TS_DiscGrad *) ts->data;
300db4653c2SDaniel Finn   PetscReal      norm, shift = 1/(0.5*ts->time_step);
301db4653c2SDaniel Finn   PetscInt       n;
302db4653c2SDaniel Finn   Vec            X0, Xdot, Xp, Xdiff;
303db4653c2SDaniel Finn   Mat            S;
304db4653c2SDaniel Finn   PetscScalar    F=0, F0=0, Gp;
305db4653c2SDaniel Finn   Vec            G, SgF;
30640be0ff1SMatthew G. Knepley   DM             dm, dmsave;
30740be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
30840be0ff1SMatthew G. Knepley 
30940be0ff1SMatthew G. Knepley   PetscFunctionBegin;
31040be0ff1SMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
31140be0ff1SMatthew G. Knepley 
312db4653c2SDaniel Finn   ierr = VecDuplicate(y, &Xp);CHKERRQ(ierr);
313db4653c2SDaniel Finn   ierr = VecDuplicate(y, &Xdiff);CHKERRQ(ierr);
314db4653c2SDaniel Finn   ierr = VecDuplicate(y, &SgF);CHKERRQ(ierr);
315db4653c2SDaniel Finn   ierr = VecDuplicate(y, &G);CHKERRQ(ierr);
316db4653c2SDaniel Finn 
317db4653c2SDaniel Finn   ierr = VecGetLocalSize(y, &n);CHKERRQ(ierr);
318db4653c2SDaniel Finn   ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
319db4653c2SDaniel Finn   ierr = MatSetSizes(S, PETSC_DECIDE, PETSC_DECIDE, n, n);CHKERRQ(ierr);
320db4653c2SDaniel Finn   ierr = MatSetFromOptions(S);CHKERRQ(ierr);
321db4653c2SDaniel Finn   ierr = MatSetUp(S);CHKERRQ(ierr);
322db4653c2SDaniel Finn   ierr = MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
323db4653c2SDaniel Finn   ierr = MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
324db4653c2SDaniel Finn 
325db4653c2SDaniel Finn   ierr = TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
326db4653c2SDaniel Finn   ierr = VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x);CHKERRQ(ierr); /* Xdot = shift (x - X0) */
327db4653c2SDaniel Finn 
328db4653c2SDaniel Finn   ierr = VecAXPBYPCZ(Xp, -1, 2, 0, X0, x);CHKERRQ(ierr); /* Xp = 2*x - X0 + (0)*Xmid */
329db4653c2SDaniel Finn   ierr = VecAXPBYPCZ(Xdiff, -1, 1, 0, X0, Xp);CHKERRQ(ierr); /* Xdiff = xp - X0 + (0)*Xdiff */
330db4653c2SDaniel Finn 
331db4653c2SDaniel Finn   if (dg->gonzalez) {
332db4653c2SDaniel Finn     ierr = (*dg->Sfunc)(ts, dg->stage_time, x ,   S,  dg->funcCtx);CHKERRQ(ierr);
333db4653c2SDaniel Finn     ierr = (*dg->Ffunc)(ts, dg->stage_time, Xp,  &F,  dg->funcCtx);CHKERRQ(ierr);
334db4653c2SDaniel Finn     ierr = (*dg->Ffunc)(ts, dg->stage_time, X0,  &F0, dg->funcCtx);CHKERRQ(ierr);
335db4653c2SDaniel Finn     ierr = (*dg->Gfunc)(ts, dg->stage_time, x ,   G,  dg->funcCtx);CHKERRQ(ierr);
336db4653c2SDaniel Finn 
337db4653c2SDaniel Finn     /* Adding Extra Gonzalez Term */
338db4653c2SDaniel Finn     ierr = VecDot(Xdiff, G, &Gp);CHKERRQ(ierr);
339db4653c2SDaniel Finn     ierr = VecNorm(Xdiff, NORM_2, &norm);CHKERRQ(ierr);
340db4653c2SDaniel Finn     if (norm < PETSC_SQRT_MACHINE_EPSILON) {
341db4653c2SDaniel Finn       Gp = 0;
342db4653c2SDaniel Finn     } else {
343db4653c2SDaniel Finn       /* Gp = (1/|xn+1 - xn|^2) * (F(xn+1) - F(xn) - Gp) */
344db4653c2SDaniel Finn       Gp = (F - F0 - Gp)/PetscSqr(norm);
345db4653c2SDaniel Finn     }
346db4653c2SDaniel Finn     ierr = VecAXPY(G, Gp, Xdiff);CHKERRQ(ierr);
347db4653c2SDaniel Finn     ierr = MatMult(S, G , SgF);CHKERRQ(ierr); /* S*gradF */
348db4653c2SDaniel Finn 
349db4653c2SDaniel Finn   } else {
350db4653c2SDaniel Finn     ierr = (*dg->Sfunc)(ts, dg->stage_time, x, S,  dg->funcCtx);CHKERRQ(ierr);
351db4653c2SDaniel Finn     ierr = (*dg->Gfunc)(ts, dg->stage_time, x, G,  dg->funcCtx);CHKERRQ(ierr);
352db4653c2SDaniel Finn 
353db4653c2SDaniel Finn     ierr = MatMult(S, G , SgF);CHKERRQ(ierr); /* Xdot = S*gradF */
354db4653c2SDaniel Finn   }
35540be0ff1SMatthew G. Knepley   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
35640be0ff1SMatthew G. Knepley   dmsave = ts->dm;
35740be0ff1SMatthew G. Knepley   ts->dm = dm;
358db4653c2SDaniel Finn   ierr = VecAXPBYPCZ(y, 1, -1, 0, Xdot, SgF);CHKERRQ(ierr);
35940be0ff1SMatthew G. Knepley   ts->dm = dmsave;
36040be0ff1SMatthew G. Knepley   ierr   = TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr);
361db4653c2SDaniel Finn 
362db4653c2SDaniel Finn   ierr = VecDestroy(&Xp);CHKERRQ(ierr);
363db4653c2SDaniel Finn   ierr = VecDestroy(&Xdiff);CHKERRQ(ierr);
364db4653c2SDaniel Finn   ierr = VecDestroy(&SgF);CHKERRQ(ierr);
365db4653c2SDaniel Finn   ierr = VecDestroy(&G);CHKERRQ(ierr);
366db4653c2SDaniel Finn   ierr = MatDestroy(&S);CHKERRQ(ierr);
367db4653c2SDaniel Finn 
36840be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
36940be0ff1SMatthew G. Knepley }
37040be0ff1SMatthew G. Knepley 
37140be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts)
37240be0ff1SMatthew G. Knepley {
37340be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg    = (TS_DiscGrad *) ts->data;
37440be0ff1SMatthew G. Knepley   PetscReal      shift = 1/(0.5*ts->time_step);
37540be0ff1SMatthew G. Knepley   Vec            Xdot;
37640be0ff1SMatthew G. Knepley   DM             dm,dmsave;
37740be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
37840be0ff1SMatthew G. Knepley 
37940be0ff1SMatthew G. Knepley   PetscFunctionBegin;
38040be0ff1SMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
38140be0ff1SMatthew G. Knepley   /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */
38240be0ff1SMatthew G. Knepley   ierr = TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr);
38340be0ff1SMatthew G. Knepley 
38440be0ff1SMatthew G. Knepley   dmsave = ts->dm;
38540be0ff1SMatthew G. Knepley   ts->dm = dm;
38640be0ff1SMatthew G. Knepley   ierr   = TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE);CHKERRQ(ierr);
38740be0ff1SMatthew G. Knepley   ts->dm = dmsave;
38840be0ff1SMatthew G. Knepley   ierr   = TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr);
38940be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
39040be0ff1SMatthew G. Knepley }
39140be0ff1SMatthew G. Knepley 
392db4653c2SDaniel Finn static PetscErrorCode TSDiscGradGetFormulation_DiscGrad(TS ts, PetscErrorCode (**Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**Gfunc)(TS, PetscReal, Vec, Vec, void *), void *ctx)
39340be0ff1SMatthew G. Knepley {
39440be0ff1SMatthew G. Knepley   TS_DiscGrad *dg = (TS_DiscGrad *) ts->data;
39540be0ff1SMatthew G. Knepley 
39640be0ff1SMatthew G. Knepley   PetscFunctionBegin;
39740be0ff1SMatthew G. Knepley   *Sfunc = dg->Sfunc;
39840be0ff1SMatthew G. Knepley   *Ffunc = dg->Ffunc;
39940be0ff1SMatthew G. Knepley   *Gfunc = dg->Gfunc;
40040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
40140be0ff1SMatthew G. Knepley }
40240be0ff1SMatthew G. Knepley 
403db4653c2SDaniel Finn static PetscErrorCode TSDiscGradSetFormulation_DiscGrad(TS ts, PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *), void *ctx)
40440be0ff1SMatthew G. Knepley {
40540be0ff1SMatthew G. Knepley   TS_DiscGrad *dg = (TS_DiscGrad *) ts->data;
40640be0ff1SMatthew G. Knepley 
40740be0ff1SMatthew G. Knepley   PetscFunctionBegin;
40840be0ff1SMatthew G. Knepley   dg->Sfunc = Sfunc;
40940be0ff1SMatthew G. Knepley   dg->Ffunc = Ffunc;
41040be0ff1SMatthew G. Knepley   dg->Gfunc = Gfunc;
411db4653c2SDaniel Finn   dg->funcCtx = ctx;
41240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
41340be0ff1SMatthew G. Knepley }
41440be0ff1SMatthew G. Knepley 
41540be0ff1SMatthew G. Knepley /*MC
41640be0ff1SMatthew G. Knepley   TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method
41740be0ff1SMatthew G. Knepley 
4182ceb51e8SPatrick Sanan   Notes:
4192ceb51e8SPatrick Sanan   This is the implicit midpoint rule, with an optional term that guarantees the discrete gradient property. This
42040be0ff1SMatthew G. Knepley   timestepper applies to systems of the form
42140be0ff1SMatthew G. Knepley $ u_t = S(u) grad F(u)
42240be0ff1SMatthew G. Knepley   where S(u) is a linear operator, and F is a functional of u.
42340be0ff1SMatthew G. Knepley 
424*f1a722f8SMatthew G. Knepley   Level: intermediate
425*f1a722f8SMatthew G. Knepley 
426db4653c2SDaniel Finn .seealso: TSCreate(), TSSetType(), TS, TSDISCGRAD, TSDiscGradSetFormulation()
42740be0ff1SMatthew G. Knepley M*/
42840be0ff1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts)
42940be0ff1SMatthew G. Knepley {
43040be0ff1SMatthew G. Knepley   TS_DiscGrad       *th;
43140be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
43240be0ff1SMatthew G. Knepley 
43340be0ff1SMatthew G. Knepley   PetscFunctionBegin;
43440be0ff1SMatthew G. Knepley   ierr = PetscCitationsRegister(DGCitation, &DGCite);CHKERRQ(ierr);
43540be0ff1SMatthew G. Knepley   ts->ops->reset           = TSReset_DiscGrad;
43640be0ff1SMatthew G. Knepley   ts->ops->destroy         = TSDestroy_DiscGrad;
43740be0ff1SMatthew G. Knepley   ts->ops->view            = TSView_DiscGrad;
43840be0ff1SMatthew G. Knepley   ts->ops->setfromoptions  = TSSetFromOptions_DiscGrad;
43940be0ff1SMatthew G. Knepley   ts->ops->setup           = TSSetUp_DiscGrad;
44040be0ff1SMatthew G. Knepley   ts->ops->step            = TSStep_DiscGrad;
44140be0ff1SMatthew G. Knepley   ts->ops->interpolate     = TSInterpolate_DiscGrad;
44240be0ff1SMatthew G. Knepley   ts->ops->getstages       = TSGetStages_DiscGrad;
44340be0ff1SMatthew G. Knepley   ts->ops->snesfunction    = SNESTSFormFunction_DiscGrad;
44440be0ff1SMatthew G. Knepley   ts->ops->snesjacobian    = SNESTSFormJacobian_DiscGrad;
44540be0ff1SMatthew G. Knepley   ts->default_adapt_type   = TSADAPTNONE;
44640be0ff1SMatthew G. Knepley 
44740be0ff1SMatthew G. Knepley   ts->usessnes = PETSC_TRUE;
44840be0ff1SMatthew G. Knepley 
44940be0ff1SMatthew G. Knepley   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
45040be0ff1SMatthew G. Knepley   ts->data = (void*)th;
451db4653c2SDaniel Finn 
452db4653c2SDaniel Finn   th->gonzalez = PETSC_FALSE;
453db4653c2SDaniel Finn 
45440be0ff1SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",TSDiscGradGetFormulation_DiscGrad);CHKERRQ(ierr);
45540be0ff1SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",TSDiscGradSetFormulation_DiscGrad);CHKERRQ(ierr);
456db4653c2SDaniel Finn   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradIsGonzalez_C",TSDiscGradIsGonzalez_DiscGrad);CHKERRQ(ierr);
457db4653c2SDaniel Finn   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradUseGonzalez_C",TSDiscGradUseGonzalez_DiscGrad);CHKERRQ(ierr);
45840be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
45940be0ff1SMatthew G. Knepley }
46040be0ff1SMatthew G. Knepley 
46140be0ff1SMatthew G. Knepley /*@C
46240be0ff1SMatthew G. Knepley   TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the formulation u_t = S grad F
46340be0ff1SMatthew G. Knepley 
46440be0ff1SMatthew G. Knepley   Not Collective
46540be0ff1SMatthew G. Knepley 
46640be0ff1SMatthew G. Knepley   Input Parameter:
46740be0ff1SMatthew G. Knepley . ts - timestepping context
46840be0ff1SMatthew G. Knepley 
46940be0ff1SMatthew G. Knepley   Output Parameters:
47040be0ff1SMatthew G. Knepley + Sfunc - constructor for the S matrix from the formulation
47140be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation
472*f1a722f8SMatthew G. Knepley . Gfunc - constructor for the gradient of F from the formulation
473*f1a722f8SMatthew G. Knepley - ctx   - the user context
47440be0ff1SMatthew G. Knepley 
47540be0ff1SMatthew G. Knepley   Calling sequence of Sfunc:
47640be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *)
47740be0ff1SMatthew G. Knepley 
47840be0ff1SMatthew G. Knepley   Calling sequence of Ffunc:
47940be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *)
48040be0ff1SMatthew G. Knepley 
48140be0ff1SMatthew G. Knepley   Calling sequence of Gfunc:
48240be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *)
48340be0ff1SMatthew G. Knepley 
48440be0ff1SMatthew G. Knepley   Level: intermediate
48540be0ff1SMatthew G. Knepley 
48640be0ff1SMatthew G. Knepley .seealso: TSDiscGradSetFormulation()
48740be0ff1SMatthew G. Knepley @*/
488db4653c2SDaniel Finn PetscErrorCode TSDiscGradGetFormulation(TS ts, PetscErrorCode (**Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**Gfunc)(TS, PetscReal, Vec, Vec, void *), void *ctx)
48940be0ff1SMatthew G. Knepley {
49040be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
49140be0ff1SMatthew G. Knepley 
49240be0ff1SMatthew G. Knepley   PetscFunctionBegin;
49340be0ff1SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
49440be0ff1SMatthew G. Knepley   PetscValidPointer(Sfunc, 2);
49540be0ff1SMatthew G. Knepley   PetscValidPointer(Ffunc, 3);
49640be0ff1SMatthew G. Knepley   PetscValidPointer(Gfunc, 4);
497db4653c2SDaniel Finn   ierr = PetscUseMethod(ts,"TSDiscGradGetFormulation_C",(TS,PetscErrorCode(**Sfunc)(TS,PetscReal,Vec,Mat,void*),PetscErrorCode(**Ffunc)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode(**Gfunc)(TS,PetscReal,Vec,Vec,void*), void*),(ts,Sfunc,Ffunc,Gfunc,ctx));CHKERRQ(ierr);
49840be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
49940be0ff1SMatthew G. Knepley }
50040be0ff1SMatthew G. Knepley 
50140be0ff1SMatthew G. Knepley /*@C
50240be0ff1SMatthew G. Knepley   TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u)
50340be0ff1SMatthew G. Knepley 
50440be0ff1SMatthew G. Knepley   Not Collective
50540be0ff1SMatthew G. Knepley 
50640be0ff1SMatthew G. Knepley   Input Parameters:
50740be0ff1SMatthew G. Knepley + ts    - timestepping context
50840be0ff1SMatthew G. Knepley . Sfunc - constructor for the S matrix from the formulation
50940be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation
51040be0ff1SMatthew G. Knepley - Gfunc - constructor for the gradient of F from the formulation
51140be0ff1SMatthew G. Knepley   Calling sequence of Sfunc:
51240be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *)
51340be0ff1SMatthew G. Knepley 
51440be0ff1SMatthew G. Knepley   Calling sequence of Ffunc:
51540be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *)
51640be0ff1SMatthew G. Knepley 
51740be0ff1SMatthew G. Knepley   Calling sequence of Gfunc:
51840be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *)
51940be0ff1SMatthew G. Knepley 
52040be0ff1SMatthew G. Knepley   Level: Intermediate
52140be0ff1SMatthew G. Knepley 
52240be0ff1SMatthew G. Knepley .seealso: TSDiscGradGetFormulation()
52340be0ff1SMatthew G. Knepley @*/
52440be0ff1SMatthew G. Knepley PetscErrorCode TSDiscGradSetFormulation(TS ts, PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*Ffunc)(TS, PetscReal, Vec , PetscScalar *, void *), PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *), void *ctx)
52540be0ff1SMatthew G. Knepley {
52640be0ff1SMatthew G. Knepley   PetscErrorCode ierr;
52740be0ff1SMatthew G. Knepley 
52840be0ff1SMatthew G. Knepley   PetscFunctionBegin;
52940be0ff1SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
53040be0ff1SMatthew G. Knepley   PetscValidFunction(Sfunc, 2);
53140be0ff1SMatthew G. Knepley   PetscValidFunction(Ffunc, 3);
53240be0ff1SMatthew G. Knepley   PetscValidFunction(Gfunc, 4);
533db4653c2SDaniel Finn   ierr = PetscTryMethod(ts,"TSDiscGradSetFormulation_C",(TS,PetscErrorCode(*Sfunc)(TS,PetscReal,Vec,Mat,void*),PetscErrorCode(*Ffunc)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode(*Gfunc)(TS,PetscReal,Vec,Vec,void*), void*),(ts,Sfunc,Ffunc,Gfunc,ctx));CHKERRQ(ierr);
534db4653c2SDaniel Finn   PetscFunctionReturn(0);
535db4653c2SDaniel Finn }
536db4653c2SDaniel Finn 
537db4653c2SDaniel Finn /*@
538db4653c2SDaniel Finn   TSDiscGradIsGonzalez - Checks flag for whether to use additional conservative terms in discrete gradient formulation.
539db4653c2SDaniel Finn 
540db4653c2SDaniel Finn   Not Collective
541db4653c2SDaniel Finn 
542db4653c2SDaniel Finn   Input Parameter:
543db4653c2SDaniel Finn .  ts - timestepping context
544db4653c2SDaniel Finn 
545db4653c2SDaniel Finn   Output Parameter:
546db4653c2SDaniel Finn .  gonzalez - PETSC_TRUE when using the Gonzalez term
547db4653c2SDaniel Finn 
548db4653c2SDaniel Finn   Level: Advanced
549db4653c2SDaniel Finn 
550db4653c2SDaniel Finn .seealso: TSDiscGradUseGonzalez(), TSDISCGRAD
551db4653c2SDaniel Finn @*/
552db4653c2SDaniel Finn PetscErrorCode TSDiscGradIsGonzalez(TS ts,PetscBool *gonzalez)
553db4653c2SDaniel Finn {
554db4653c2SDaniel Finn   PetscErrorCode ierr;
555db4653c2SDaniel Finn 
556db4653c2SDaniel Finn   PetscFunctionBegin;
557db4653c2SDaniel Finn   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
558db4653c2SDaniel Finn   PetscValidPointer(gonzalez,2);
559db4653c2SDaniel Finn   ierr = PetscUseMethod(ts,"TSDiscGradIsGonzalez_C",(TS,PetscBool*),(ts,gonzalez));CHKERRQ(ierr);
560db4653c2SDaniel Finn   PetscFunctionReturn(0);
561db4653c2SDaniel Finn }
562db4653c2SDaniel Finn 
563db4653c2SDaniel Finn /*@
564db4653c2SDaniel Finn   TSDiscGradUseGonzalez - Sets discrete gradient formulation with or without additional conservative terms.  Without flag, the discrete gradients timestepper is just backwards euler
565db4653c2SDaniel Finn 
566db4653c2SDaniel Finn   Not Collective
567db4653c2SDaniel Finn 
568*f1a722f8SMatthew G. Knepley   Input Parameters:
569db4653c2SDaniel Finn + ts - timestepping context
570db4653c2SDaniel Finn - flg - PETSC_TRUE to use the Gonzalez term
571db4653c2SDaniel Finn 
572db4653c2SDaniel Finn   Options Database:
573db4653c2SDaniel Finn . -ts_discgrad_gonzalez <flg>
574db4653c2SDaniel Finn 
575db4653c2SDaniel Finn   Level: Intermediate
576db4653c2SDaniel Finn 
577db4653c2SDaniel Finn .seealso: TSDISCGRAD
578db4653c2SDaniel Finn @*/
579db4653c2SDaniel Finn PetscErrorCode TSDiscGradUseGonzalez(TS ts,PetscBool flg)
580db4653c2SDaniel Finn {
581db4653c2SDaniel Finn   PetscErrorCode ierr;
582db4653c2SDaniel Finn 
583db4653c2SDaniel Finn   PetscFunctionBegin;
584db4653c2SDaniel Finn   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
585db4653c2SDaniel Finn   ierr = PetscTryMethod(ts,"TSDiscGradUseGonzalez_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
58640be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
58740be0ff1SMatthew G. Knepley }
588