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; 2040be0ff1SMatthew G. Knepley PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *); 2140be0ff1SMatthew G. Knepley PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *); 2240be0ff1SMatthew G. Knepley PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *); 2340be0ff1SMatthew G. Knepley } TS_DiscGrad; 2440be0ff1SMatthew G. Knepley 2540be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGradGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 2640be0ff1SMatthew G. Knepley { 2740be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 2840be0ff1SMatthew G. Knepley PetscErrorCode ierr; 2940be0ff1SMatthew G. Knepley 3040be0ff1SMatthew G. Knepley PetscFunctionBegin; 3140be0ff1SMatthew G. Knepley if (X0) { 3240be0ff1SMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0);CHKERRQ(ierr);} 3340be0ff1SMatthew G. Knepley else {*X0 = ts->vec_sol;} 3440be0ff1SMatthew G. Knepley } 3540be0ff1SMatthew G. Knepley if (Xdot) { 3640be0ff1SMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot);CHKERRQ(ierr);} 3740be0ff1SMatthew G. Knepley else {*Xdot = dg->Xdot;} 3840be0ff1SMatthew G. Knepley } 3940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 4040be0ff1SMatthew G. Knepley } 4140be0ff1SMatthew G. Knepley 4240be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 4340be0ff1SMatthew G. Knepley { 4440be0ff1SMatthew G. Knepley PetscErrorCode ierr; 4540be0ff1SMatthew G. Knepley 4640be0ff1SMatthew G. Knepley PetscFunctionBegin; 4740be0ff1SMatthew G. Knepley if (X0) { 4840be0ff1SMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0);CHKERRQ(ierr);} 4940be0ff1SMatthew G. Knepley } 5040be0ff1SMatthew G. Knepley if (Xdot) { 5140be0ff1SMatthew G. Knepley if (dm && dm != ts->dm) {ierr = DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot);CHKERRQ(ierr);} 5240be0ff1SMatthew G. Knepley } 5340be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 5440be0ff1SMatthew G. Knepley } 5540be0ff1SMatthew G. Knepley 5640be0ff1SMatthew G. Knepley static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, void *ctx) 5740be0ff1SMatthew G. Knepley { 5840be0ff1SMatthew G. Knepley PetscFunctionBegin; 5940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 6040be0ff1SMatthew G. Knepley } 6140be0ff1SMatthew G. Knepley 6240be0ff1SMatthew G. Knepley static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 6340be0ff1SMatthew G. Knepley { 6440be0ff1SMatthew G. Knepley TS ts = (TS) ctx; 6540be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_c, Xdot_c; 6640be0ff1SMatthew G. Knepley PetscErrorCode ierr; 6740be0ff1SMatthew G. Knepley 6840be0ff1SMatthew G. Knepley PetscFunctionBegin; 6940be0ff1SMatthew G. Knepley ierr = TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot);CHKERRQ(ierr); 7040be0ff1SMatthew G. Knepley ierr = TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c);CHKERRQ(ierr); 7140be0ff1SMatthew G. Knepley ierr = MatRestrict(restrct, X0, X0_c);CHKERRQ(ierr); 7240be0ff1SMatthew G. Knepley ierr = MatRestrict(restrct, Xdot, Xdot_c);CHKERRQ(ierr); 7340be0ff1SMatthew G. Knepley ierr = VecPointwiseMult(X0_c, rscale, X0_c);CHKERRQ(ierr); 7440be0ff1SMatthew G. Knepley ierr = VecPointwiseMult(Xdot_c, rscale, Xdot_c);CHKERRQ(ierr); 7540be0ff1SMatthew G. Knepley ierr = TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot);CHKERRQ(ierr); 7640be0ff1SMatthew G. Knepley ierr = TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c);CHKERRQ(ierr); 7740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 7840be0ff1SMatthew G. Knepley } 7940be0ff1SMatthew G. Knepley 8040be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, void *ctx) 8140be0ff1SMatthew G. Knepley { 8240be0ff1SMatthew G. Knepley PetscFunctionBegin; 8340be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 8440be0ff1SMatthew G. Knepley } 8540be0ff1SMatthew G. Knepley 8640be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 8740be0ff1SMatthew G. Knepley { 8840be0ff1SMatthew G. Knepley TS ts = (TS) ctx; 8940be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_sub, Xdot_sub; 9040be0ff1SMatthew G. Knepley PetscErrorCode ierr; 9140be0ff1SMatthew G. Knepley 9240be0ff1SMatthew G. Knepley PetscFunctionBegin; 9340be0ff1SMatthew G. Knepley ierr = TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 9440be0ff1SMatthew G. Knepley ierr = TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub);CHKERRQ(ierr); 9540be0ff1SMatthew G. Knepley 9640be0ff1SMatthew G. Knepley ierr = VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 9740be0ff1SMatthew G. Knepley ierr = VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 9840be0ff1SMatthew G. Knepley 9940be0ff1SMatthew G. Knepley ierr = VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 10040be0ff1SMatthew G. Knepley ierr = VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 10140be0ff1SMatthew G. Knepley 10240be0ff1SMatthew G. Knepley ierr = TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 10340be0ff1SMatthew G. Knepley ierr = TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub);CHKERRQ(ierr); 10440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 10540be0ff1SMatthew G. Knepley } 10640be0ff1SMatthew G. Knepley 10740be0ff1SMatthew G. Knepley static PetscErrorCode TSSetUp_DiscGrad(TS ts) 10840be0ff1SMatthew G. Knepley { 10940be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 11040be0ff1SMatthew G. Knepley DM dm; 11140be0ff1SMatthew G. Knepley PetscErrorCode ierr; 11240be0ff1SMatthew G. Knepley 11340be0ff1SMatthew G. Knepley PetscFunctionBegin; 11440be0ff1SMatthew G. Knepley if (!dg->X) {ierr = VecDuplicate(ts->vec_sol, &dg->X);CHKERRQ(ierr);} 11540be0ff1SMatthew G. Knepley if (!dg->X0) {ierr = VecDuplicate(ts->vec_sol, &dg->X0);CHKERRQ(ierr);} 11640be0ff1SMatthew G. Knepley if (!dg->Xdot) {ierr = VecDuplicate(ts->vec_sol, &dg->Xdot);CHKERRQ(ierr);} 11740be0ff1SMatthew G. Knepley 11840be0ff1SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 11940be0ff1SMatthew G. Knepley ierr = DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr); 12040be0ff1SMatthew G. Knepley ierr = DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr); 12140be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 12240be0ff1SMatthew G. Knepley } 12340be0ff1SMatthew G. Knepley 12440be0ff1SMatthew G. Knepley static PetscErrorCode TSSetFromOptions_DiscGrad(PetscOptionItems *PetscOptionsObject, TS ts) 12540be0ff1SMatthew G. Knepley { 12640be0ff1SMatthew G. Knepley PetscErrorCode ierr; 12740be0ff1SMatthew G. Knepley 12840be0ff1SMatthew G. Knepley PetscFunctionBegin; 12940be0ff1SMatthew G. Knepley ierr = PetscOptionsHead(PetscOptionsObject, "Discrete Gradients ODE solver options");CHKERRQ(ierr); 13040be0ff1SMatthew G. Knepley { 13140be0ff1SMatthew G. Knepley //ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSDiscGradSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 13240be0ff1SMatthew G. Knepley } 13340be0ff1SMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 13440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 13540be0ff1SMatthew G. Knepley } 13640be0ff1SMatthew G. Knepley 13740be0ff1SMatthew G. Knepley static PetscErrorCode TSView_DiscGrad(TS ts,PetscViewer viewer) 13840be0ff1SMatthew G. Knepley { 13940be0ff1SMatthew G. Knepley PetscBool iascii; 14040be0ff1SMatthew G. Knepley PetscErrorCode ierr; 14140be0ff1SMatthew G. Knepley 14240be0ff1SMatthew G. Knepley PetscFunctionBegin; 14340be0ff1SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 14440be0ff1SMatthew G. Knepley if (iascii) { 14540be0ff1SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer," Discrete Gradients\n");CHKERRQ(ierr); 14640be0ff1SMatthew G. Knepley } 14740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 14840be0ff1SMatthew G. Knepley } 14940be0ff1SMatthew G. Knepley 15040be0ff1SMatthew G. Knepley static PetscErrorCode TSReset_DiscGrad(TS ts) 15140be0ff1SMatthew G. Knepley { 15240be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 15340be0ff1SMatthew G. Knepley PetscErrorCode ierr; 15440be0ff1SMatthew G. Knepley 15540be0ff1SMatthew G. Knepley PetscFunctionBegin; 15640be0ff1SMatthew G. Knepley ierr = VecDestroy(&dg->X);CHKERRQ(ierr); 15740be0ff1SMatthew G. Knepley ierr = VecDestroy(&dg->X0);CHKERRQ(ierr); 15840be0ff1SMatthew G. Knepley ierr = VecDestroy(&dg->Xdot);CHKERRQ(ierr); 15940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 16040be0ff1SMatthew G. Knepley } 16140be0ff1SMatthew G. Knepley 16240be0ff1SMatthew G. Knepley static PetscErrorCode TSDestroy_DiscGrad(TS ts) 16340be0ff1SMatthew G. Knepley { 16440be0ff1SMatthew G. Knepley DM dm; 16540be0ff1SMatthew G. Knepley PetscErrorCode ierr; 16640be0ff1SMatthew G. Knepley 16740be0ff1SMatthew G. Knepley PetscFunctionBegin; 16840be0ff1SMatthew G. Knepley ierr = TSReset_DiscGrad(ts);CHKERRQ(ierr); 16940be0ff1SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 17040be0ff1SMatthew G. Knepley if (dm) { 17140be0ff1SMatthew G. Knepley ierr = DMCoarsenHookRemove(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr); 17240be0ff1SMatthew G. Knepley ierr = DMSubDomainHookRemove(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts);CHKERRQ(ierr); 17340be0ff1SMatthew G. Knepley } 17440be0ff1SMatthew G. Knepley ierr = PetscFree(ts->data);CHKERRQ(ierr); 17540be0ff1SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",NULL);CHKERRQ(ierr); 17640be0ff1SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",NULL);CHKERRQ(ierr); 17740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 17840be0ff1SMatthew G. Knepley } 17940be0ff1SMatthew G. Knepley 18040be0ff1SMatthew G. Knepley static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X) 18140be0ff1SMatthew G. Knepley { 18240be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad*)ts->data; 18340be0ff1SMatthew G. Knepley PetscReal dt = t - ts->ptime; 18440be0ff1SMatthew G. Knepley PetscErrorCode ierr; 18540be0ff1SMatthew G. Knepley 18640be0ff1SMatthew G. Knepley PetscFunctionBegin; 18740be0ff1SMatthew G. Knepley ierr = VecCopy(ts->vec_sol, dg->X);CHKERRQ(ierr); 18840be0ff1SMatthew G. Knepley ierr = VecWAXPY(X, dt, dg->Xdot, dg->X);CHKERRQ(ierr); 18940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 19040be0ff1SMatthew G. Knepley } 19140be0ff1SMatthew G. Knepley 19240be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x) 19340be0ff1SMatthew G. Knepley { 19440be0ff1SMatthew G. Knepley SNES snes; 19540be0ff1SMatthew G. Knepley PetscInt nits, lits; 19640be0ff1SMatthew G. Knepley PetscErrorCode ierr; 19740be0ff1SMatthew G. Knepley 19840be0ff1SMatthew G. Knepley PetscFunctionBegin; 19940be0ff1SMatthew G. Knepley ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 20040be0ff1SMatthew G. Knepley ierr = SNESSolve(snes, b, x);CHKERRQ(ierr); 20140be0ff1SMatthew G. Knepley ierr = SNESGetIterationNumber(snes, &nits);CHKERRQ(ierr); 20240be0ff1SMatthew G. Knepley ierr = SNESGetLinearSolveIterations(snes, &lits);CHKERRQ(ierr); 20340be0ff1SMatthew G. Knepley ts->snes_its += nits; 20440be0ff1SMatthew G. Knepley ts->ksp_its += lits; 20540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 20640be0ff1SMatthew G. Knepley } 20740be0ff1SMatthew G. Knepley 20840be0ff1SMatthew G. Knepley static PetscErrorCode TSStep_DiscGrad(TS ts) 20940be0ff1SMatthew G. Knepley { 21040be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 21140be0ff1SMatthew G. Knepley TSAdapt adapt; 21240be0ff1SMatthew G. Knepley TSStepStatus status = TS_STEP_INCOMPLETE; 21340be0ff1SMatthew G. Knepley PetscInt rejections = 0; 21440be0ff1SMatthew G. Knepley PetscBool stageok, accept = PETSC_TRUE; 21540be0ff1SMatthew G. Knepley PetscReal next_time_step = ts->time_step; 21640be0ff1SMatthew G. Knepley PetscErrorCode ierr; 21740be0ff1SMatthew G. Knepley 21840be0ff1SMatthew G. Knepley PetscFunctionBegin; 21940be0ff1SMatthew G. Knepley ierr = TSGetAdapt(ts, &adapt);CHKERRQ(ierr); 22040be0ff1SMatthew G. Knepley if (!ts->steprollback) {ierr = VecCopy(ts->vec_sol, dg->X0);CHKERRQ(ierr);} 22140be0ff1SMatthew G. Knepley 22240be0ff1SMatthew G. Knepley while (!ts->reason && status != TS_STEP_COMPLETE) { 22340be0ff1SMatthew G. Knepley PetscReal shift = 1/(0.5*ts->time_step); 22440be0ff1SMatthew G. Knepley 22540be0ff1SMatthew G. Knepley dg->stage_time = ts->ptime + 0.5*ts->time_step; 22640be0ff1SMatthew G. Knepley 22740be0ff1SMatthew G. Knepley ierr = VecCopy(dg->X0, dg->X);CHKERRQ(ierr); 22840be0ff1SMatthew G. Knepley ierr = TSPreStage(ts, dg->stage_time);CHKERRQ(ierr); 22940be0ff1SMatthew G. Knepley ierr = TSDiscGrad_SNESSolve(ts, NULL, dg->X);CHKERRQ(ierr); 23040be0ff1SMatthew G. Knepley ierr = TSPostStage(ts, dg->stage_time, 0, &dg->X);CHKERRQ(ierr); 23140be0ff1SMatthew G. Knepley ierr = TSAdaptCheckStage(adapt, ts, dg->stage_time, dg->X, &stageok);CHKERRQ(ierr); 23240be0ff1SMatthew G. Knepley if (!stageok) goto reject_step; 23340be0ff1SMatthew G. Knepley 23440be0ff1SMatthew G. Knepley status = TS_STEP_PENDING; 23540be0ff1SMatthew G. Knepley ierr = VecAXPBYPCZ(dg->Xdot, -shift, shift, 0, dg->X0, dg->X);CHKERRQ(ierr); 23640be0ff1SMatthew G. Knepley ierr = VecAXPY(ts->vec_sol, ts->time_step, dg->Xdot);CHKERRQ(ierr); 23740be0ff1SMatthew G. Knepley ierr = TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept);CHKERRQ(ierr); 23840be0ff1SMatthew G. Knepley status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 23940be0ff1SMatthew G. Knepley if (!accept) { 24040be0ff1SMatthew G. Knepley ierr = VecCopy(dg->X0, ts->vec_sol);CHKERRQ(ierr); 24140be0ff1SMatthew G. Knepley ts->time_step = next_time_step; 24240be0ff1SMatthew G. Knepley goto reject_step; 24340be0ff1SMatthew G. Knepley } 24440be0ff1SMatthew G. Knepley ts->ptime += ts->time_step; 24540be0ff1SMatthew G. Knepley ts->time_step = next_time_step; 24640be0ff1SMatthew G. Knepley break; 24740be0ff1SMatthew G. Knepley 24840be0ff1SMatthew G. Knepley reject_step: 24940be0ff1SMatthew G. Knepley ts->reject++; accept = PETSC_FALSE; 25040be0ff1SMatthew G. Knepley if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) { 25140be0ff1SMatthew G. Knepley ts->reason = TS_DIVERGED_STEP_REJECTED; 25240be0ff1SMatthew G. Knepley ierr = PetscInfo2(ts, "Step=%D, step rejections %D greater than current TS allowed, stopping solve\n", ts->steps, rejections);CHKERRQ(ierr); 25340be0ff1SMatthew G. Knepley } 25440be0ff1SMatthew G. Knepley } 25540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 25640be0ff1SMatthew G. Knepley } 25740be0ff1SMatthew G. Knepley 25840be0ff1SMatthew G. Knepley static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y) 25940be0ff1SMatthew G. Knepley { 26040be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 26140be0ff1SMatthew G. Knepley 26240be0ff1SMatthew G. Knepley PetscFunctionBegin; 26340be0ff1SMatthew G. Knepley if (ns) *ns = 1; 26440be0ff1SMatthew G. Knepley if (Y) *Y = &(dg->X); 26540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 26640be0ff1SMatthew G. Knepley } 26740be0ff1SMatthew G. Knepley 26840be0ff1SMatthew G. Knepley /* 26940be0ff1SMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 27040be0ff1SMatthew G. Knepley G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0 27140be0ff1SMatthew G. Knepley */ 27240be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts) 27340be0ff1SMatthew G. Knepley { 27440be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 27540be0ff1SMatthew G. Knepley PetscReal shift = 1/(0.5*ts->time_step); 27640be0ff1SMatthew G. Knepley Vec X0, Xdot; 27740be0ff1SMatthew G. Knepley DM dm, dmsave; 27840be0ff1SMatthew G. Knepley PetscErrorCode ierr; 27940be0ff1SMatthew G. Knepley 28040be0ff1SMatthew G. Knepley PetscFunctionBegin; 28140be0ff1SMatthew G. Knepley ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 28240be0ff1SMatthew G. Knepley ierr = TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 28340be0ff1SMatthew G. Knepley ierr = VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x);CHKERRQ(ierr); 28440be0ff1SMatthew G. Knepley 28540be0ff1SMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 28640be0ff1SMatthew G. Knepley dmsave = ts->dm; 28740be0ff1SMatthew G. Knepley ts->dm = dm; 28840be0ff1SMatthew G. Knepley ierr = TSComputeIFunction(ts, dg->stage_time, x, Xdot, y, PETSC_FALSE);CHKERRQ(ierr); 28940be0ff1SMatthew G. Knepley ts->dm = dmsave; 29040be0ff1SMatthew G. Knepley ierr = TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot);CHKERRQ(ierr); 29140be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 29240be0ff1SMatthew G. Knepley } 29340be0ff1SMatthew G. Knepley 29440be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts) 29540be0ff1SMatthew G. Knepley { 29640be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 29740be0ff1SMatthew G. Knepley PetscReal shift = 1/(0.5*ts->time_step); 29840be0ff1SMatthew G. Knepley Vec Xdot; 29940be0ff1SMatthew G. Knepley DM dm,dmsave; 30040be0ff1SMatthew G. Knepley PetscErrorCode ierr; 30140be0ff1SMatthew G. Knepley 30240be0ff1SMatthew G. Knepley PetscFunctionBegin; 30340be0ff1SMatthew G. Knepley ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 30440be0ff1SMatthew G. Knepley /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */ 30540be0ff1SMatthew G. Knepley ierr = TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr); 30640be0ff1SMatthew G. Knepley 30740be0ff1SMatthew G. Knepley dmsave = ts->dm; 30840be0ff1SMatthew G. Knepley ts->dm = dm; 30940be0ff1SMatthew G. Knepley ierr = TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE);CHKERRQ(ierr); 31040be0ff1SMatthew G. Knepley ts->dm = dmsave; 31140be0ff1SMatthew G. Knepley ierr = TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot);CHKERRQ(ierr); 31240be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 31340be0ff1SMatthew G. Knepley } 31440be0ff1SMatthew G. Knepley 31540be0ff1SMatthew G. Knepley 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 *)) 31640be0ff1SMatthew G. Knepley { 31740be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 31840be0ff1SMatthew G. Knepley 31940be0ff1SMatthew G. Knepley PetscFunctionBegin; 32040be0ff1SMatthew G. Knepley *Sfunc = dg->Sfunc; 32140be0ff1SMatthew G. Knepley *Ffunc = dg->Ffunc; 32240be0ff1SMatthew G. Knepley *Gfunc = dg->Gfunc; 32340be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 32440be0ff1SMatthew G. Knepley } 32540be0ff1SMatthew G. Knepley 32640be0ff1SMatthew G. Knepley 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 *)) 32740be0ff1SMatthew G. Knepley { 32840be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 32940be0ff1SMatthew G. Knepley 33040be0ff1SMatthew G. Knepley PetscFunctionBegin; 33140be0ff1SMatthew G. Knepley dg->Sfunc = Sfunc; 33240be0ff1SMatthew G. Knepley dg->Ffunc = Ffunc; 33340be0ff1SMatthew G. Knepley dg->Gfunc = Gfunc; 33440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 33540be0ff1SMatthew G. Knepley } 33640be0ff1SMatthew G. Knepley 33740be0ff1SMatthew G. Knepley /*MC 33840be0ff1SMatthew G. Knepley TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method 33940be0ff1SMatthew G. Knepley 340*2ceb51e8SPatrick Sanan Notes: 341*2ceb51e8SPatrick Sanan This is the implicit midpoint rule, with an optional term that guarantees the discrete gradient property. This 34240be0ff1SMatthew G. Knepley timestepper applies to systems of the form 34340be0ff1SMatthew G. Knepley $ u_t = S(u) grad F(u) 34440be0ff1SMatthew G. Knepley where S(u) is a linear operator, and F is a functional of u. 34540be0ff1SMatthew G. Knepley 346*2ceb51e8SPatrick Sanan Reference: 347*2ceb51e8SPatrick Sanan Gonzalez, Oscar. "Time integration and discrete Hamiltonian systems." Journal of Nonlinear Science 6.5 (1996): 449-467. 348*2ceb51e8SPatrick Sanan 349*2ceb51e8SPatrick Sanan Level: beginner 350*2ceb51e8SPatrick Sanan 35140be0ff1SMatthew G. Knepley .seealso: TSCreate(), TSSetType(), TS, TSTHETA, TSDiscGradSetFormulation() 35240be0ff1SMatthew G. Knepley M*/ 35340be0ff1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts) 35440be0ff1SMatthew G. Knepley { 35540be0ff1SMatthew G. Knepley TS_DiscGrad *th; 35640be0ff1SMatthew G. Knepley PetscErrorCode ierr; 35740be0ff1SMatthew G. Knepley 35840be0ff1SMatthew G. Knepley PetscFunctionBegin; 35940be0ff1SMatthew G. Knepley ierr = PetscCitationsRegister(DGCitation, &DGCite);CHKERRQ(ierr); 36040be0ff1SMatthew G. Knepley ts->ops->reset = TSReset_DiscGrad; 36140be0ff1SMatthew G. Knepley ts->ops->destroy = TSDestroy_DiscGrad; 36240be0ff1SMatthew G. Knepley ts->ops->view = TSView_DiscGrad; 36340be0ff1SMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_DiscGrad; 36440be0ff1SMatthew G. Knepley ts->ops->setup = TSSetUp_DiscGrad; 36540be0ff1SMatthew G. Knepley ts->ops->step = TSStep_DiscGrad; 36640be0ff1SMatthew G. Knepley ts->ops->interpolate = TSInterpolate_DiscGrad; 36740be0ff1SMatthew G. Knepley ts->ops->getstages = TSGetStages_DiscGrad; 36840be0ff1SMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_DiscGrad; 36940be0ff1SMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_DiscGrad; 37040be0ff1SMatthew G. Knepley ts->default_adapt_type = TSADAPTNONE; 37140be0ff1SMatthew G. Knepley 37240be0ff1SMatthew G. Knepley ts->usessnes = PETSC_TRUE; 37340be0ff1SMatthew G. Knepley 37440be0ff1SMatthew G. Knepley ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 37540be0ff1SMatthew G. Knepley ts->data = (void*)th; 37640be0ff1SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",TSDiscGradGetFormulation_DiscGrad);CHKERRQ(ierr); 37740be0ff1SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",TSDiscGradSetFormulation_DiscGrad);CHKERRQ(ierr); 37840be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 37940be0ff1SMatthew G. Knepley } 38040be0ff1SMatthew G. Knepley 38140be0ff1SMatthew G. Knepley /*@C 38240be0ff1SMatthew G. Knepley TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the formulation u_t = S grad F 38340be0ff1SMatthew G. Knepley 38440be0ff1SMatthew G. Knepley Not Collective 38540be0ff1SMatthew G. Knepley 38640be0ff1SMatthew G. Knepley Input Parameter: 38740be0ff1SMatthew G. Knepley . ts - timestepping context 38840be0ff1SMatthew G. Knepley 38940be0ff1SMatthew G. Knepley Output Parameters: 39040be0ff1SMatthew G. Knepley + Sfunc - constructor for the S matrix from the formulation 39140be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 39240be0ff1SMatthew G. Knepley - Gfunc - constructor for the gradient of F from the formulation 39340be0ff1SMatthew G. Knepley 39440be0ff1SMatthew G. Knepley Calling sequence of Sfunc: 39540be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 39640be0ff1SMatthew G. Knepley 39740be0ff1SMatthew G. Knepley Calling sequence of Ffunc: 39840be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 39940be0ff1SMatthew G. Knepley 40040be0ff1SMatthew G. Knepley Calling sequence of Gfunc: 40140be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 40240be0ff1SMatthew G. Knepley 40340be0ff1SMatthew G. Knepley Level: intermediate 40440be0ff1SMatthew G. Knepley 40540be0ff1SMatthew G. Knepley .seealso: TSDiscGradSetFormulation() 40640be0ff1SMatthew G. Knepley @*/ 40740be0ff1SMatthew G. Knepley 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 *)) 40840be0ff1SMatthew G. Knepley { 40940be0ff1SMatthew G. Knepley PetscErrorCode ierr; 41040be0ff1SMatthew G. Knepley 41140be0ff1SMatthew G. Knepley PetscFunctionBegin; 41240be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 41340be0ff1SMatthew G. Knepley PetscValidPointer(Sfunc, 2); 41440be0ff1SMatthew G. Knepley PetscValidPointer(Ffunc, 3); 41540be0ff1SMatthew G. Knepley PetscValidPointer(Gfunc, 4); 41640be0ff1SMatthew G. Knepley 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*)),(ts,Sfunc,Ffunc,Gfunc));CHKERRQ(ierr); 41740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 41840be0ff1SMatthew G. Knepley } 41940be0ff1SMatthew G. Knepley 42040be0ff1SMatthew G. Knepley /*@C 42140be0ff1SMatthew G. Knepley TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u) 42240be0ff1SMatthew G. Knepley 42340be0ff1SMatthew G. Knepley Not Collective 42440be0ff1SMatthew G. Knepley 42540be0ff1SMatthew G. Knepley Input Parameters: 42640be0ff1SMatthew G. Knepley + ts - timestepping context 42740be0ff1SMatthew G. Knepley . Sfunc - constructor for the S matrix from the formulation 42840be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 42940be0ff1SMatthew G. Knepley - Gfunc - constructor for the gradient of F from the formulation 43040be0ff1SMatthew G. Knepley 43140be0ff1SMatthew G. Knepley Calling sequence of Sfunc: 43240be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 43340be0ff1SMatthew G. Knepley 43440be0ff1SMatthew G. Knepley Calling sequence of Ffunc: 43540be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 43640be0ff1SMatthew G. Knepley 43740be0ff1SMatthew G. Knepley Calling sequence of Gfunc: 43840be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 43940be0ff1SMatthew G. Knepley 44040be0ff1SMatthew G. Knepley Level: Intermediate 44140be0ff1SMatthew G. Knepley 44240be0ff1SMatthew G. Knepley .seealso: TSDiscGradGetFormulation() 44340be0ff1SMatthew G. Knepley @*/ 44440be0ff1SMatthew 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) 44540be0ff1SMatthew G. Knepley { 44640be0ff1SMatthew G. Knepley PetscErrorCode ierr; 44740be0ff1SMatthew G. Knepley 44840be0ff1SMatthew G. Knepley PetscFunctionBegin; 44940be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 45040be0ff1SMatthew G. Knepley PetscValidFunction(Sfunc, 2); 45140be0ff1SMatthew G. Knepley PetscValidFunction(Ffunc, 3); 45240be0ff1SMatthew G. Knepley PetscValidFunction(Gfunc, 4); 45340be0ff1SMatthew G. Knepley 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*)),(ts,Sfunc,Ffunc,Gfunc));CHKERRQ(ierr); 45440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 45540be0ff1SMatthew G. Knepley } 456