xref: /petsc/src/ts/impls/implicit/discgrad/tsdiscgrad.c (revision dadcf80911fb48939c55327431ae8d7e47dbe367)
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 
3140be0ff1SMatthew G. Knepley   PetscFunctionBegin;
3240be0ff1SMatthew G. Knepley   if (X0) {
335f80ce2aSJacob Faibussowitsch     if (dm && dm != ts->dm) CHKERRQ(DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0));
3440be0ff1SMatthew G. Knepley     else                    {*X0  = ts->vec_sol;}
3540be0ff1SMatthew G. Knepley   }
3640be0ff1SMatthew G. Knepley   if (Xdot) {
375f80ce2aSJacob Faibussowitsch     if (dm && dm != ts->dm) CHKERRQ(DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot));
3840be0ff1SMatthew G. Knepley     else                    {*Xdot = dg->Xdot;}
3940be0ff1SMatthew G. Knepley   }
4040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
4140be0ff1SMatthew G. Knepley }
4240be0ff1SMatthew G. Knepley 
4340be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
4440be0ff1SMatthew G. Knepley {
4540be0ff1SMatthew G. Knepley   PetscFunctionBegin;
4640be0ff1SMatthew G. Knepley   if (X0) {
475f80ce2aSJacob Faibussowitsch     if (dm && dm != ts->dm) CHKERRQ(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0));
4840be0ff1SMatthew G. Knepley   }
4940be0ff1SMatthew G. Knepley   if (Xdot) {
505f80ce2aSJacob Faibussowitsch     if (dm && dm != ts->dm) CHKERRQ(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot));
5140be0ff1SMatthew G. Knepley   }
5240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
5340be0ff1SMatthew G. Knepley }
5440be0ff1SMatthew G. Knepley 
5540be0ff1SMatthew G. Knepley static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, void *ctx)
5640be0ff1SMatthew G. Knepley {
5740be0ff1SMatthew G. Knepley   PetscFunctionBegin;
5840be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
5940be0ff1SMatthew G. Knepley }
6040be0ff1SMatthew G. Knepley 
6140be0ff1SMatthew G. Knepley static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
6240be0ff1SMatthew G. Knepley {
6340be0ff1SMatthew G. Knepley   TS             ts = (TS) ctx;
6440be0ff1SMatthew G. Knepley   Vec            X0, Xdot, X0_c, Xdot_c;
6540be0ff1SMatthew G. Knepley 
6640be0ff1SMatthew G. Knepley   PetscFunctionBegin;
675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRestrict(restrct, X0, X0_c));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRestrict(restrct, Xdot, Xdot_c));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(X0_c, rscale, X0_c));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(Xdot_c, rscale, Xdot_c));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
7540be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
7640be0ff1SMatthew G. Knepley }
7740be0ff1SMatthew G. Knepley 
7840be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, void *ctx)
7940be0ff1SMatthew G. Knepley {
8040be0ff1SMatthew G. Knepley   PetscFunctionBegin;
8140be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
8240be0ff1SMatthew G. Knepley }
8340be0ff1SMatthew G. Knepley 
8440be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
8540be0ff1SMatthew G. Knepley {
8640be0ff1SMatthew G. Knepley   TS             ts = (TS) ctx;
8740be0ff1SMatthew G. Knepley   Vec            X0, Xdot, X0_sub, Xdot_sub;
8840be0ff1SMatthew G. Knepley 
8940be0ff1SMatthew G. Knepley   PetscFunctionBegin;
905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
9240be0ff1SMatthew G. Knepley 
935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
9540be0ff1SMatthew G. Knepley 
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
9840be0ff1SMatthew G. Knepley 
995f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
10140be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
10240be0ff1SMatthew G. Knepley }
10340be0ff1SMatthew G. Knepley 
10440be0ff1SMatthew G. Knepley static PetscErrorCode TSSetUp_DiscGrad(TS ts)
10540be0ff1SMatthew G. Knepley {
10640be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
10740be0ff1SMatthew G. Knepley   DM             dm;
10840be0ff1SMatthew G. Knepley 
10940be0ff1SMatthew G. Knepley   PetscFunctionBegin;
1105f80ce2aSJacob Faibussowitsch   if (!dg->X)    CHKERRQ(VecDuplicate(ts->vec_sol, &dg->X));
1115f80ce2aSJacob Faibussowitsch   if (!dg->X0)   CHKERRQ(VecDuplicate(ts->vec_sol, &dg->X0));
1125f80ce2aSJacob Faibussowitsch   if (!dg->Xdot) CHKERRQ(VecDuplicate(ts->vec_sol, &dg->Xdot));
11340be0ff1SMatthew G. Knepley 
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts));
11740be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
11840be0ff1SMatthew G. Knepley }
11940be0ff1SMatthew G. Knepley 
12040be0ff1SMatthew G. Knepley static PetscErrorCode TSSetFromOptions_DiscGrad(PetscOptionItems *PetscOptionsObject, TS ts)
12140be0ff1SMatthew G. Knepley {
122db4653c2SDaniel Finn   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
12340be0ff1SMatthew G. Knepley 
12440be0ff1SMatthew G. Knepley   PetscFunctionBegin;
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHead(PetscOptionsObject, "Discrete Gradients ODE solver options"));
12640be0ff1SMatthew G. Knepley   {
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsBool("-ts_discgrad_gonzalez","Use Gonzalez term in discrete gradients formulation","TSDiscGradUseGonzalez",dg->gonzalez,&dg->gonzalez,NULL));
12840be0ff1SMatthew G. Knepley   }
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsTail());
13040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
13140be0ff1SMatthew G. Knepley }
13240be0ff1SMatthew G. Knepley 
13340be0ff1SMatthew G. Knepley static PetscErrorCode TSView_DiscGrad(TS ts,PetscViewer viewer)
13440be0ff1SMatthew G. Knepley {
13540be0ff1SMatthew G. Knepley   PetscBool      iascii;
13640be0ff1SMatthew G. Knepley 
13740be0ff1SMatthew G. Knepley   PetscFunctionBegin;
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
13940be0ff1SMatthew G. Knepley   if (iascii) {
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Discrete Gradients\n"));
14140be0ff1SMatthew G. Knepley   }
14240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
14340be0ff1SMatthew G. Knepley }
14440be0ff1SMatthew G. Knepley 
145db4653c2SDaniel Finn static PetscErrorCode TSDiscGradIsGonzalez_DiscGrad(TS ts,PetscBool *gonzalez)
146db4653c2SDaniel Finn {
147db4653c2SDaniel Finn   TS_DiscGrad *dg = (TS_DiscGrad*)ts->data;
148db4653c2SDaniel Finn 
149db4653c2SDaniel Finn   PetscFunctionBegin;
150db4653c2SDaniel Finn   *gonzalez = dg->gonzalez;
151db4653c2SDaniel Finn   PetscFunctionReturn(0);
152db4653c2SDaniel Finn }
153db4653c2SDaniel Finn 
154db4653c2SDaniel Finn static PetscErrorCode TSDiscGradUseGonzalez_DiscGrad(TS ts,PetscBool flg)
155db4653c2SDaniel Finn {
156db4653c2SDaniel Finn   TS_DiscGrad *dg = (TS_DiscGrad*)ts->data;
157db4653c2SDaniel Finn 
158db4653c2SDaniel Finn   PetscFunctionBegin;
159db4653c2SDaniel Finn   dg->gonzalez = flg;
160db4653c2SDaniel Finn   PetscFunctionReturn(0);
161db4653c2SDaniel Finn }
162db4653c2SDaniel Finn 
16340be0ff1SMatthew G. Knepley static PetscErrorCode TSReset_DiscGrad(TS ts)
16440be0ff1SMatthew G. Knepley {
16540be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
16640be0ff1SMatthew G. Knepley 
16740be0ff1SMatthew G. Knepley   PetscFunctionBegin;
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&dg->X));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&dg->X0));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&dg->Xdot));
17140be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
17240be0ff1SMatthew G. Knepley }
17340be0ff1SMatthew G. Knepley 
17440be0ff1SMatthew G. Knepley static PetscErrorCode TSDestroy_DiscGrad(TS ts)
17540be0ff1SMatthew G. Knepley {
17640be0ff1SMatthew G. Knepley   DM             dm;
17740be0ff1SMatthew G. Knepley 
17840be0ff1SMatthew G. Knepley   PetscFunctionBegin;
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(TSReset_DiscGrad(ts));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
18140be0ff1SMatthew G. Knepley   if (dm) {
1825f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCoarsenHookRemove(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts));
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSubDomainHookRemove(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts));
18440be0ff1SMatthew G. Knepley   }
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ts->data));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",NULL));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",NULL));
18840be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
18940be0ff1SMatthew G. Knepley }
19040be0ff1SMatthew G. Knepley 
19140be0ff1SMatthew G. Knepley static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X)
19240be0ff1SMatthew G. Knepley {
19340be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad*)ts->data;
19440be0ff1SMatthew G. Knepley   PetscReal      dt = t - ts->ptime;
19540be0ff1SMatthew G. Knepley 
19640be0ff1SMatthew G. Knepley   PetscFunctionBegin;
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(ts->vec_sol, dg->X));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(X, dt, dg->Xdot, dg->X));
19940be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
20040be0ff1SMatthew G. Knepley }
20140be0ff1SMatthew G. Knepley 
20240be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x)
20340be0ff1SMatthew G. Knepley {
20440be0ff1SMatthew G. Knepley   SNES           snes;
20540be0ff1SMatthew G. Knepley   PetscInt       nits, lits;
20640be0ff1SMatthew G. Knepley 
20740be0ff1SMatthew G. Knepley   PetscFunctionBegin;
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSNES(ts, &snes));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes, b, x));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes, &nits));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLinearSolveIterations(snes, &lits));
21240be0ff1SMatthew G. Knepley   ts->snes_its += nits;
21340be0ff1SMatthew G. Knepley   ts->ksp_its  += lits;
21440be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
21540be0ff1SMatthew G. Knepley }
21640be0ff1SMatthew G. Knepley 
21740be0ff1SMatthew G. Knepley static PetscErrorCode TSStep_DiscGrad(TS ts)
21840be0ff1SMatthew G. Knepley {
21940be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg = (TS_DiscGrad *) ts->data;
22040be0ff1SMatthew G. Knepley   TSAdapt        adapt;
22140be0ff1SMatthew G. Knepley   TSStepStatus   status          = TS_STEP_INCOMPLETE;
22240be0ff1SMatthew G. Knepley   PetscInt       rejections      = 0;
22340be0ff1SMatthew G. Knepley   PetscBool      stageok, accept = PETSC_TRUE;
22440be0ff1SMatthew G. Knepley   PetscReal      next_time_step  = ts->time_step;
22540be0ff1SMatthew G. Knepley 
22640be0ff1SMatthew G. Knepley   PetscFunctionBegin;
2275f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetAdapt(ts, &adapt));
2285f80ce2aSJacob Faibussowitsch   if (!ts->steprollback) CHKERRQ(VecCopy(ts->vec_sol, dg->X0));
22940be0ff1SMatthew G. Knepley 
23040be0ff1SMatthew G. Knepley   while (!ts->reason && status != TS_STEP_COMPLETE) {
23140be0ff1SMatthew G. Knepley     PetscReal shift = 1/(0.5*ts->time_step);
23240be0ff1SMatthew G. Knepley 
23340be0ff1SMatthew G. Knepley     dg->stage_time = ts->ptime + 0.5*ts->time_step;
23440be0ff1SMatthew G. Knepley 
2355f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(dg->X0, dg->X));
2365f80ce2aSJacob Faibussowitsch     CHKERRQ(TSPreStage(ts, dg->stage_time));
2375f80ce2aSJacob Faibussowitsch     CHKERRQ(TSDiscGrad_SNESSolve(ts, NULL, dg->X));
2385f80ce2aSJacob Faibussowitsch     CHKERRQ(TSPostStage(ts, dg->stage_time, 0, &dg->X));
2395f80ce2aSJacob Faibussowitsch     CHKERRQ(TSAdaptCheckStage(adapt, ts, dg->stage_time, dg->X, &stageok));
24040be0ff1SMatthew G. Knepley     if (!stageok) goto reject_step;
24140be0ff1SMatthew G. Knepley 
24240be0ff1SMatthew G. Knepley     status = TS_STEP_PENDING;
2435f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPBYPCZ(dg->Xdot, -shift, shift, 0, dg->X0, dg->X));
2445f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(ts->vec_sol, ts->time_step, dg->Xdot));
2455f80ce2aSJacob Faibussowitsch     CHKERRQ(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
24640be0ff1SMatthew G. Knepley     status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
24740be0ff1SMatthew G. Knepley     if (!accept) {
2485f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(dg->X0, ts->vec_sol));
24940be0ff1SMatthew G. Knepley       ts->time_step = next_time_step;
25040be0ff1SMatthew G. Knepley       goto reject_step;
25140be0ff1SMatthew G. Knepley     }
25240be0ff1SMatthew G. Knepley     ts->ptime    += ts->time_step;
25340be0ff1SMatthew G. Knepley     ts->time_step = next_time_step;
25440be0ff1SMatthew G. Knepley     break;
25540be0ff1SMatthew G. Knepley 
25640be0ff1SMatthew G. Knepley   reject_step:
25740be0ff1SMatthew G. Knepley     ts->reject++; accept = PETSC_FALSE;
25840be0ff1SMatthew G. Knepley     if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) {
25940be0ff1SMatthew G. Knepley       ts->reason = TS_DIVERGED_STEP_REJECTED;
2605f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(ts, "Step=%D, step rejections %D greater than current TS allowed, stopping solve\n", ts->steps, rejections));
26140be0ff1SMatthew G. Knepley     }
26240be0ff1SMatthew G. Knepley   }
26340be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
26440be0ff1SMatthew G. Knepley }
26540be0ff1SMatthew G. Knepley 
26640be0ff1SMatthew G. Knepley static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y)
26740be0ff1SMatthew G. Knepley {
26840be0ff1SMatthew G. Knepley   TS_DiscGrad *dg = (TS_DiscGrad *) ts->data;
26940be0ff1SMatthew G. Knepley 
27040be0ff1SMatthew G. Knepley   PetscFunctionBegin;
27140be0ff1SMatthew G. Knepley   if (ns) *ns = 1;
27240be0ff1SMatthew G. Knepley   if (Y)  *Y  = &(dg->X);
27340be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
27440be0ff1SMatthew G. Knepley }
27540be0ff1SMatthew G. Knepley 
27640be0ff1SMatthew G. Knepley /*
27740be0ff1SMatthew G. Knepley   This defines the nonlinear equation that is to be solved with SNES
27840be0ff1SMatthew G. Knepley     G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0
27940be0ff1SMatthew G. Knepley */
280db4653c2SDaniel Finn 
281db4653c2SDaniel Finn /* x = (x+x')/2 */
282db4653c2SDaniel Finn /* NEED TO CALCULATE x_{n+1} from x and x_{n}*/
28340be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts)
28440be0ff1SMatthew G. Knepley {
285db4653c2SDaniel Finn 
28640be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg    = (TS_DiscGrad *) ts->data;
287db4653c2SDaniel Finn   PetscReal      norm, shift = 1/(0.5*ts->time_step);
288db4653c2SDaniel Finn   PetscInt       n;
289db4653c2SDaniel Finn   Vec            X0, Xdot, Xp, Xdiff;
290db4653c2SDaniel Finn   Mat            S;
291db4653c2SDaniel Finn   PetscScalar    F=0, F0=0, Gp;
292db4653c2SDaniel Finn   Vec            G, SgF;
29340be0ff1SMatthew G. Knepley   DM             dm, dmsave;
29440be0ff1SMatthew G. Knepley 
29540be0ff1SMatthew G. Knepley   PetscFunctionBegin;
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes, &dm));
29740be0ff1SMatthew G. Knepley 
2985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(y, &Xp));
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(y, &Xdiff));
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(y, &SgF));
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(y, &G));
302db4653c2SDaniel Finn 
3035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(y, &n));
3045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &S));
3055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(S, PETSC_DECIDE, PETSC_DECIDE, n, n));
3065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(S));
3075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(S));
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY));
310db4653c2SDaniel Finn 
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot));
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); /* Xdot = shift (x - X0) */
313db4653c2SDaniel Finn 
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPBYPCZ(Xp, -1, 2, 0, X0, x)); /* Xp = 2*x - X0 + (0)*Xmid */
3155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPBYPCZ(Xdiff, -1, 1, 0, X0, Xp)); /* Xdiff = xp - X0 + (0)*Xdiff */
316db4653c2SDaniel Finn 
317db4653c2SDaniel Finn   if (dg->gonzalez) {
3185f80ce2aSJacob Faibussowitsch     CHKERRQ((*dg->Sfunc)(ts, dg->stage_time, x ,   S,  dg->funcCtx));
3195f80ce2aSJacob Faibussowitsch     CHKERRQ((*dg->Ffunc)(ts, dg->stage_time, Xp,  &F,  dg->funcCtx));
3205f80ce2aSJacob Faibussowitsch     CHKERRQ((*dg->Ffunc)(ts, dg->stage_time, X0,  &F0, dg->funcCtx));
3215f80ce2aSJacob Faibussowitsch     CHKERRQ((*dg->Gfunc)(ts, dg->stage_time, x ,   G,  dg->funcCtx));
322db4653c2SDaniel Finn 
323db4653c2SDaniel Finn     /* Adding Extra Gonzalez Term */
3245f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDot(Xdiff, G, &Gp));
3255f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(Xdiff, NORM_2, &norm));
326db4653c2SDaniel Finn     if (norm < PETSC_SQRT_MACHINE_EPSILON) {
327db4653c2SDaniel Finn       Gp = 0;
328db4653c2SDaniel Finn     } else {
329db4653c2SDaniel Finn       /* Gp = (1/|xn+1 - xn|^2) * (F(xn+1) - F(xn) - Gp) */
330db4653c2SDaniel Finn       Gp = (F - F0 - Gp)/PetscSqr(norm);
331db4653c2SDaniel Finn     }
3325f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(G, Gp, Xdiff));
3335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(S, G , SgF)); /* S*gradF */
334db4653c2SDaniel Finn 
335db4653c2SDaniel Finn   } else {
3365f80ce2aSJacob Faibussowitsch     CHKERRQ((*dg->Sfunc)(ts, dg->stage_time, x, S,  dg->funcCtx));
3375f80ce2aSJacob Faibussowitsch     CHKERRQ((*dg->Gfunc)(ts, dg->stage_time, x, G,  dg->funcCtx));
338db4653c2SDaniel Finn 
3395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(S, G , SgF)); /* Xdot = S*gradF */
340db4653c2SDaniel Finn   }
34140be0ff1SMatthew G. Knepley   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
34240be0ff1SMatthew G. Knepley   dmsave = ts->dm;
34340be0ff1SMatthew G. Knepley   ts->dm = dm;
3445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPBYPCZ(y, 1, -1, 0, Xdot, SgF));
34540be0ff1SMatthew G. Knepley   ts->dm = dmsave;
3465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot));
347db4653c2SDaniel Finn 
3485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Xp));
3495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Xdiff));
3505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&SgF));
3515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&G));
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&S));
353db4653c2SDaniel Finn 
35440be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
35540be0ff1SMatthew G. Knepley }
35640be0ff1SMatthew G. Knepley 
35740be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts)
35840be0ff1SMatthew G. Knepley {
35940be0ff1SMatthew G. Knepley   TS_DiscGrad   *dg    = (TS_DiscGrad *) ts->data;
36040be0ff1SMatthew G. Knepley   PetscReal      shift = 1/(0.5*ts->time_step);
36140be0ff1SMatthew G. Knepley   Vec            Xdot;
36240be0ff1SMatthew G. Knepley   DM             dm,dmsave;
36340be0ff1SMatthew G. Knepley 
36440be0ff1SMatthew G. Knepley   PetscFunctionBegin;
3655f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes, &dm));
36640be0ff1SMatthew G. Knepley   /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot));
36840be0ff1SMatthew G. Knepley 
36940be0ff1SMatthew G. Knepley   dmsave = ts->dm;
37040be0ff1SMatthew G. Knepley   ts->dm = dm;
3715f80ce2aSJacob Faibussowitsch   CHKERRQ(TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
37240be0ff1SMatthew G. Knepley   ts->dm = dmsave;
3735f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot));
37440be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
37540be0ff1SMatthew G. Knepley }
37640be0ff1SMatthew G. Knepley 
377db4653c2SDaniel 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)
37840be0ff1SMatthew G. Knepley {
37940be0ff1SMatthew G. Knepley   TS_DiscGrad *dg = (TS_DiscGrad *) ts->data;
38040be0ff1SMatthew G. Knepley 
38140be0ff1SMatthew G. Knepley   PetscFunctionBegin;
38240be0ff1SMatthew G. Knepley   *Sfunc = dg->Sfunc;
38340be0ff1SMatthew G. Knepley   *Ffunc = dg->Ffunc;
38440be0ff1SMatthew G. Knepley   *Gfunc = dg->Gfunc;
38540be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
38640be0ff1SMatthew G. Knepley }
38740be0ff1SMatthew G. Knepley 
388db4653c2SDaniel 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)
38940be0ff1SMatthew G. Knepley {
39040be0ff1SMatthew G. Knepley   TS_DiscGrad *dg = (TS_DiscGrad *) ts->data;
39140be0ff1SMatthew G. Knepley 
39240be0ff1SMatthew G. Knepley   PetscFunctionBegin;
39340be0ff1SMatthew G. Knepley   dg->Sfunc = Sfunc;
39440be0ff1SMatthew G. Knepley   dg->Ffunc = Ffunc;
39540be0ff1SMatthew G. Knepley   dg->Gfunc = Gfunc;
396db4653c2SDaniel Finn   dg->funcCtx = ctx;
39740be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
39840be0ff1SMatthew G. Knepley }
39940be0ff1SMatthew G. Knepley 
40040be0ff1SMatthew G. Knepley /*MC
40140be0ff1SMatthew G. Knepley   TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method
40240be0ff1SMatthew G. Knepley 
4032ceb51e8SPatrick Sanan   Notes:
4042ceb51e8SPatrick Sanan   This is the implicit midpoint rule, with an optional term that guarantees the discrete gradient property. This
40540be0ff1SMatthew G. Knepley   timestepper applies to systems of the form
40640be0ff1SMatthew G. Knepley $ u_t = S(u) grad F(u)
40740be0ff1SMatthew G. Knepley   where S(u) is a linear operator, and F is a functional of u.
40840be0ff1SMatthew G. Knepley 
409f1a722f8SMatthew G. Knepley   Level: intermediate
410f1a722f8SMatthew G. Knepley 
411db4653c2SDaniel Finn .seealso: TSCreate(), TSSetType(), TS, TSDISCGRAD, TSDiscGradSetFormulation()
41240be0ff1SMatthew G. Knepley M*/
41340be0ff1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts)
41440be0ff1SMatthew G. Knepley {
41540be0ff1SMatthew G. Knepley   TS_DiscGrad       *th;
41640be0ff1SMatthew G. Knepley 
41740be0ff1SMatthew G. Knepley   PetscFunctionBegin;
4185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCitationsRegister(DGCitation, &DGCite));
41940be0ff1SMatthew G. Knepley   ts->ops->reset           = TSReset_DiscGrad;
42040be0ff1SMatthew G. Knepley   ts->ops->destroy         = TSDestroy_DiscGrad;
42140be0ff1SMatthew G. Knepley   ts->ops->view            = TSView_DiscGrad;
42240be0ff1SMatthew G. Knepley   ts->ops->setfromoptions  = TSSetFromOptions_DiscGrad;
42340be0ff1SMatthew G. Knepley   ts->ops->setup           = TSSetUp_DiscGrad;
42440be0ff1SMatthew G. Knepley   ts->ops->step            = TSStep_DiscGrad;
42540be0ff1SMatthew G. Knepley   ts->ops->interpolate     = TSInterpolate_DiscGrad;
42640be0ff1SMatthew G. Knepley   ts->ops->getstages       = TSGetStages_DiscGrad;
42740be0ff1SMatthew G. Knepley   ts->ops->snesfunction    = SNESTSFormFunction_DiscGrad;
42840be0ff1SMatthew G. Knepley   ts->ops->snesjacobian    = SNESTSFormJacobian_DiscGrad;
42940be0ff1SMatthew G. Knepley   ts->default_adapt_type   = TSADAPTNONE;
43040be0ff1SMatthew G. Knepley 
43140be0ff1SMatthew G. Knepley   ts->usessnes = PETSC_TRUE;
43240be0ff1SMatthew G. Knepley 
4335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(ts,&th));
43440be0ff1SMatthew G. Knepley   ts->data = (void*)th;
435db4653c2SDaniel Finn 
436db4653c2SDaniel Finn   th->gonzalez = PETSC_FALSE;
437db4653c2SDaniel Finn 
4385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",TSDiscGradGetFormulation_DiscGrad));
4395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",TSDiscGradSetFormulation_DiscGrad));
4405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradIsGonzalez_C",TSDiscGradIsGonzalez_DiscGrad));
4415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradUseGonzalez_C",TSDiscGradUseGonzalez_DiscGrad));
44240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
44340be0ff1SMatthew G. Knepley }
44440be0ff1SMatthew G. Knepley 
44540be0ff1SMatthew G. Knepley /*@C
44640be0ff1SMatthew G. Knepley   TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the formulation u_t = S grad F
44740be0ff1SMatthew G. Knepley 
44840be0ff1SMatthew G. Knepley   Not Collective
44940be0ff1SMatthew G. Knepley 
45040be0ff1SMatthew G. Knepley   Input Parameter:
45140be0ff1SMatthew G. Knepley . ts - timestepping context
45240be0ff1SMatthew G. Knepley 
45340be0ff1SMatthew G. Knepley   Output Parameters:
45440be0ff1SMatthew G. Knepley + Sfunc - constructor for the S matrix from the formulation
45540be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation
456f1a722f8SMatthew G. Knepley . Gfunc - constructor for the gradient of F from the formulation
457f1a722f8SMatthew G. Knepley - ctx   - the user context
45840be0ff1SMatthew G. Knepley 
45940be0ff1SMatthew G. Knepley   Calling sequence of Sfunc:
46040be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *)
46140be0ff1SMatthew G. Knepley 
46240be0ff1SMatthew G. Knepley   Calling sequence of Ffunc:
46340be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *)
46440be0ff1SMatthew G. Knepley 
46540be0ff1SMatthew G. Knepley   Calling sequence of Gfunc:
46640be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *)
46740be0ff1SMatthew G. Knepley 
46840be0ff1SMatthew G. Knepley   Level: intermediate
46940be0ff1SMatthew G. Knepley 
47040be0ff1SMatthew G. Knepley .seealso: TSDiscGradSetFormulation()
47140be0ff1SMatthew G. Knepley @*/
472db4653c2SDaniel 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)
47340be0ff1SMatthew G. Knepley {
47440be0ff1SMatthew G. Knepley   PetscFunctionBegin;
47540be0ff1SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
47640be0ff1SMatthew G. Knepley   PetscValidPointer(Sfunc, 2);
47740be0ff1SMatthew G. Knepley   PetscValidPointer(Ffunc, 3);
47840be0ff1SMatthew G. Knepley   PetscValidPointer(Gfunc, 4);
4795f80ce2aSJacob Faibussowitsch   CHKERRQ(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)));
48040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
48140be0ff1SMatthew G. Knepley }
48240be0ff1SMatthew G. Knepley 
48340be0ff1SMatthew G. Knepley /*@C
48440be0ff1SMatthew G. Knepley   TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u)
48540be0ff1SMatthew G. Knepley 
48640be0ff1SMatthew G. Knepley   Not Collective
48740be0ff1SMatthew G. Knepley 
48840be0ff1SMatthew G. Knepley   Input Parameters:
48940be0ff1SMatthew G. Knepley + ts    - timestepping context
49040be0ff1SMatthew G. Knepley . Sfunc - constructor for the S matrix from the formulation
49140be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation
49240be0ff1SMatthew G. Knepley - Gfunc - constructor for the gradient of F from the formulation
49340be0ff1SMatthew G. Knepley   Calling sequence of Sfunc:
49440be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *)
49540be0ff1SMatthew G. Knepley 
49640be0ff1SMatthew G. Knepley   Calling sequence of Ffunc:
49740be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *)
49840be0ff1SMatthew G. Knepley 
49940be0ff1SMatthew G. Knepley   Calling sequence of Gfunc:
50040be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *)
50140be0ff1SMatthew G. Knepley 
50240be0ff1SMatthew G. Knepley   Level: Intermediate
50340be0ff1SMatthew G. Knepley 
50440be0ff1SMatthew G. Knepley .seealso: TSDiscGradGetFormulation()
50540be0ff1SMatthew G. Knepley @*/
50640be0ff1SMatthew 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)
50740be0ff1SMatthew G. Knepley {
50840be0ff1SMatthew G. Knepley   PetscFunctionBegin;
50940be0ff1SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
51040be0ff1SMatthew G. Knepley   PetscValidFunction(Sfunc, 2);
51140be0ff1SMatthew G. Knepley   PetscValidFunction(Ffunc, 3);
51240be0ff1SMatthew G. Knepley   PetscValidFunction(Gfunc, 4);
5135f80ce2aSJacob Faibussowitsch   CHKERRQ(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)));
514db4653c2SDaniel Finn   PetscFunctionReturn(0);
515db4653c2SDaniel Finn }
516db4653c2SDaniel Finn 
517db4653c2SDaniel Finn /*@
518db4653c2SDaniel Finn   TSDiscGradIsGonzalez - Checks flag for whether to use additional conservative terms in discrete gradient formulation.
519db4653c2SDaniel Finn 
520db4653c2SDaniel Finn   Not Collective
521db4653c2SDaniel Finn 
522db4653c2SDaniel Finn   Input Parameter:
523db4653c2SDaniel Finn .  ts - timestepping context
524db4653c2SDaniel Finn 
525db4653c2SDaniel Finn   Output Parameter:
526db4653c2SDaniel Finn .  gonzalez - PETSC_TRUE when using the Gonzalez term
527db4653c2SDaniel Finn 
528db4653c2SDaniel Finn   Level: Advanced
529db4653c2SDaniel Finn 
530db4653c2SDaniel Finn .seealso: TSDiscGradUseGonzalez(), TSDISCGRAD
531db4653c2SDaniel Finn @*/
532db4653c2SDaniel Finn PetscErrorCode TSDiscGradIsGonzalez(TS ts,PetscBool *gonzalez)
533db4653c2SDaniel Finn {
534db4653c2SDaniel Finn   PetscFunctionBegin;
535db4653c2SDaniel Finn   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
536*dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(gonzalez,2);
5375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscUseMethod(ts,"TSDiscGradIsGonzalez_C",(TS,PetscBool*),(ts,gonzalez)));
538db4653c2SDaniel Finn   PetscFunctionReturn(0);
539db4653c2SDaniel Finn }
540db4653c2SDaniel Finn 
541db4653c2SDaniel Finn /*@
542db4653c2SDaniel Finn   TSDiscGradUseGonzalez - Sets discrete gradient formulation with or without additional conservative terms.  Without flag, the discrete gradients timestepper is just backwards euler
543db4653c2SDaniel Finn 
544db4653c2SDaniel Finn   Not Collective
545db4653c2SDaniel Finn 
546f1a722f8SMatthew G. Knepley   Input Parameters:
547db4653c2SDaniel Finn + ts - timestepping context
548db4653c2SDaniel Finn - flg - PETSC_TRUE to use the Gonzalez term
549db4653c2SDaniel Finn 
550db4653c2SDaniel Finn   Options Database:
55167b8a455SSatish Balay . -ts_discgrad_gonzalez <flg> - use the Gonzalez term for the discrete gradient formulation
552db4653c2SDaniel Finn 
553db4653c2SDaniel Finn   Level: Intermediate
554db4653c2SDaniel Finn 
555db4653c2SDaniel Finn .seealso: TSDISCGRAD
556db4653c2SDaniel Finn @*/
557db4653c2SDaniel Finn PetscErrorCode TSDiscGradUseGonzalez(TS ts,PetscBool flg)
558db4653c2SDaniel Finn {
559db4653c2SDaniel Finn   PetscFunctionBegin;
560db4653c2SDaniel Finn   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(ts,"TSDiscGradUseGonzalez_C",(TS,PetscBool),(ts,flg)));
56240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
56340be0ff1SMatthew G. Knepley }
564