1 /* 2 Code for timestepping with discrete gradient integrators 3 */ 4 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5 #include <petscdm.h> 6 7 PetscBool DGCite = PETSC_FALSE; 8 const char DGCitation[] = "@article{Gonzalez1996,\n" 9 " title = {Time integration and discrete Hamiltonian systems},\n" 10 " author = {Oscar Gonzalez},\n" 11 " journal = {Journal of Nonlinear Science},\n" 12 " volume = {6},\n" 13 " pages = {449--467},\n" 14 " doi = {10.1007/978-1-4612-1246-1_10},\n" 15 " year = {1996}\n}\n"; 16 17 typedef struct { 18 PetscReal stage_time; 19 Vec X0, X, Xdot; 20 void *funcCtx; 21 PetscBool gonzalez; 22 PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *); 23 PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *); 24 PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *); 25 } TS_DiscGrad; 26 27 static PetscErrorCode TSDiscGradGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 28 { 29 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 30 31 PetscFunctionBegin; 32 if (X0) { 33 if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 34 else {*X0 = ts->vec_sol;} 35 } 36 if (Xdot) { 37 if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 38 else {*Xdot = dg->Xdot;} 39 } 40 PetscFunctionReturn(0); 41 } 42 43 static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 44 { 45 PetscFunctionBegin; 46 if (X0) { 47 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 48 } 49 if (Xdot) { 50 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 51 } 52 PetscFunctionReturn(0); 53 } 54 55 static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, void *ctx) 56 { 57 PetscFunctionBegin; 58 PetscFunctionReturn(0); 59 } 60 61 static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 62 { 63 TS ts = (TS) ctx; 64 Vec X0, Xdot, X0_c, Xdot_c; 65 66 PetscFunctionBegin; 67 PetscCall(TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot)); 68 PetscCall(TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 69 PetscCall(MatRestrict(restrct, X0, X0_c)); 70 PetscCall(MatRestrict(restrct, Xdot, Xdot_c)); 71 PetscCall(VecPointwiseMult(X0_c, rscale, X0_c)); 72 PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c)); 73 PetscCall(TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot)); 74 PetscCall(TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 75 PetscFunctionReturn(0); 76 } 77 78 static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, void *ctx) 79 { 80 PetscFunctionBegin; 81 PetscFunctionReturn(0); 82 } 83 84 static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 85 { 86 TS ts = (TS) ctx; 87 Vec X0, Xdot, X0_sub, Xdot_sub; 88 89 PetscFunctionBegin; 90 PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 91 PetscCall(TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 92 93 PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 94 PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 95 96 PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 97 PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 98 99 PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 100 PetscCall(TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 101 PetscFunctionReturn(0); 102 } 103 104 static PetscErrorCode TSSetUp_DiscGrad(TS ts) 105 { 106 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 107 DM dm; 108 109 PetscFunctionBegin; 110 if (!dg->X) PetscCall(VecDuplicate(ts->vec_sol, &dg->X)); 111 if (!dg->X0) PetscCall(VecDuplicate(ts->vec_sol, &dg->X0)); 112 if (!dg->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &dg->Xdot)); 113 114 PetscCall(TSGetDM(ts, &dm)); 115 PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts)); 116 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts)); 117 PetscFunctionReturn(0); 118 } 119 120 static PetscErrorCode TSSetFromOptions_DiscGrad(PetscOptionItems *PetscOptionsObject, TS ts) 121 { 122 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 123 124 PetscFunctionBegin; 125 PetscOptionsHeadBegin(PetscOptionsObject, "Discrete Gradients ODE solver options"); 126 { 127 PetscCall(PetscOptionsBool("-ts_discgrad_gonzalez","Use Gonzalez term in discrete gradients formulation","TSDiscGradUseGonzalez",dg->gonzalez,&dg->gonzalez,NULL)); 128 } 129 PetscOptionsHeadEnd(); 130 PetscFunctionReturn(0); 131 } 132 133 static PetscErrorCode TSView_DiscGrad(TS ts,PetscViewer viewer) 134 { 135 PetscBool iascii; 136 137 PetscFunctionBegin; 138 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 139 if (iascii) { 140 PetscCall(PetscViewerASCIIPrintf(viewer," Discrete Gradients\n")); 141 } 142 PetscFunctionReturn(0); 143 } 144 145 static PetscErrorCode TSDiscGradIsGonzalez_DiscGrad(TS ts,PetscBool *gonzalez) 146 { 147 TS_DiscGrad *dg = (TS_DiscGrad*)ts->data; 148 149 PetscFunctionBegin; 150 *gonzalez = dg->gonzalez; 151 PetscFunctionReturn(0); 152 } 153 154 static PetscErrorCode TSDiscGradUseGonzalez_DiscGrad(TS ts,PetscBool flg) 155 { 156 TS_DiscGrad *dg = (TS_DiscGrad*)ts->data; 157 158 PetscFunctionBegin; 159 dg->gonzalez = flg; 160 PetscFunctionReturn(0); 161 } 162 163 static PetscErrorCode TSReset_DiscGrad(TS ts) 164 { 165 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 166 167 PetscFunctionBegin; 168 PetscCall(VecDestroy(&dg->X)); 169 PetscCall(VecDestroy(&dg->X0)); 170 PetscCall(VecDestroy(&dg->Xdot)); 171 PetscFunctionReturn(0); 172 } 173 174 static PetscErrorCode TSDestroy_DiscGrad(TS ts) 175 { 176 DM dm; 177 178 PetscFunctionBegin; 179 PetscCall(TSReset_DiscGrad(ts)); 180 PetscCall(TSGetDM(ts, &dm)); 181 if (dm) { 182 PetscCall(DMCoarsenHookRemove(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts)); 183 PetscCall(DMSubDomainHookRemove(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts)); 184 } 185 PetscCall(PetscFree(ts->data)); 186 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",NULL)); 187 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",NULL)); 188 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradIsGonzalez_C",NULL)); 189 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradUseGonzalez_C",NULL)); 190 PetscFunctionReturn(0); 191 } 192 193 static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X) 194 { 195 TS_DiscGrad *dg = (TS_DiscGrad*)ts->data; 196 PetscReal dt = t - ts->ptime; 197 198 PetscFunctionBegin; 199 PetscCall(VecCopy(ts->vec_sol, dg->X)); 200 PetscCall(VecWAXPY(X, dt, dg->Xdot, dg->X)); 201 PetscFunctionReturn(0); 202 } 203 204 static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x) 205 { 206 SNES snes; 207 PetscInt nits, lits; 208 209 PetscFunctionBegin; 210 PetscCall(TSGetSNES(ts, &snes)); 211 PetscCall(SNESSolve(snes, b, x)); 212 PetscCall(SNESGetIterationNumber(snes, &nits)); 213 PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 214 ts->snes_its += nits; 215 ts->ksp_its += lits; 216 PetscFunctionReturn(0); 217 } 218 219 static PetscErrorCode TSStep_DiscGrad(TS ts) 220 { 221 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 222 TSAdapt adapt; 223 TSStepStatus status = TS_STEP_INCOMPLETE; 224 PetscInt rejections = 0; 225 PetscBool stageok, accept = PETSC_TRUE; 226 PetscReal next_time_step = ts->time_step; 227 228 PetscFunctionBegin; 229 PetscCall(TSGetAdapt(ts, &adapt)); 230 if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, dg->X0)); 231 232 while (!ts->reason && status != TS_STEP_COMPLETE) { 233 PetscReal shift = 1/(0.5*ts->time_step); 234 235 dg->stage_time = ts->ptime + 0.5*ts->time_step; 236 237 PetscCall(VecCopy(dg->X0, dg->X)); 238 PetscCall(TSPreStage(ts, dg->stage_time)); 239 PetscCall(TSDiscGrad_SNESSolve(ts, NULL, dg->X)); 240 PetscCall(TSPostStage(ts, dg->stage_time, 0, &dg->X)); 241 PetscCall(TSAdaptCheckStage(adapt, ts, dg->stage_time, dg->X, &stageok)); 242 if (!stageok) goto reject_step; 243 244 status = TS_STEP_PENDING; 245 PetscCall(VecAXPBYPCZ(dg->Xdot, -shift, shift, 0, dg->X0, dg->X)); 246 PetscCall(VecAXPY(ts->vec_sol, ts->time_step, dg->Xdot)); 247 PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 248 status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 249 if (!accept) { 250 PetscCall(VecCopy(dg->X0, ts->vec_sol)); 251 ts->time_step = next_time_step; 252 goto reject_step; 253 } 254 ts->ptime += ts->time_step; 255 ts->time_step = next_time_step; 256 break; 257 258 reject_step: 259 ts->reject++; accept = PETSC_FALSE; 260 if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) { 261 ts->reason = TS_DIVERGED_STEP_REJECTED; 262 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 263 } 264 } 265 PetscFunctionReturn(0); 266 } 267 268 static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y) 269 { 270 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 271 272 PetscFunctionBegin; 273 if (ns) *ns = 1; 274 if (Y) *Y = &(dg->X); 275 PetscFunctionReturn(0); 276 } 277 278 /* 279 This defines the nonlinear equation that is to be solved with SNES 280 G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0 281 */ 282 283 /* x = (x+x')/2 */ 284 /* NEED TO CALCULATE x_{n+1} from x and x_{n}*/ 285 static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts) 286 { 287 288 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 289 PetscReal norm, shift = 1/(0.5*ts->time_step); 290 PetscInt n; 291 Vec X0, Xdot, Xp, Xdiff; 292 Mat S; 293 PetscScalar F=0, F0=0, Gp; 294 Vec G, SgF; 295 DM dm, dmsave; 296 297 PetscFunctionBegin; 298 PetscCall(SNESGetDM(snes, &dm)); 299 300 PetscCall(VecDuplicate(y, &Xp)); 301 PetscCall(VecDuplicate(y, &Xdiff)); 302 PetscCall(VecDuplicate(y, &SgF)); 303 PetscCall(VecDuplicate(y, &G)); 304 305 PetscCall(VecGetLocalSize(y, &n)); 306 PetscCall(MatCreate(PETSC_COMM_WORLD, &S)); 307 PetscCall(MatSetSizes(S, PETSC_DECIDE, PETSC_DECIDE, n, n)); 308 PetscCall(MatSetFromOptions(S)); 309 PetscCall(MatSetUp(S)); 310 PetscCall(MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY)); 311 PetscCall(MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY)); 312 313 PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 314 PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); /* Xdot = shift (x - X0) */ 315 316 PetscCall(VecAXPBYPCZ(Xp, -1, 2, 0, X0, x)); /* Xp = 2*x - X0 + (0)*Xmid */ 317 PetscCall(VecAXPBYPCZ(Xdiff, -1, 1, 0, X0, Xp)); /* Xdiff = xp - X0 + (0)*Xdiff */ 318 319 if (dg->gonzalez) { 320 PetscCall((*dg->Sfunc)(ts, dg->stage_time, x , S, dg->funcCtx)); 321 PetscCall((*dg->Ffunc)(ts, dg->stage_time, Xp, &F, dg->funcCtx)); 322 PetscCall((*dg->Ffunc)(ts, dg->stage_time, X0, &F0, dg->funcCtx)); 323 PetscCall((*dg->Gfunc)(ts, dg->stage_time, x , G, dg->funcCtx)); 324 325 /* Adding Extra Gonzalez Term */ 326 PetscCall(VecDot(Xdiff, G, &Gp)); 327 PetscCall(VecNorm(Xdiff, NORM_2, &norm)); 328 if (norm < PETSC_SQRT_MACHINE_EPSILON) { 329 Gp = 0; 330 } else { 331 /* Gp = (1/|xn+1 - xn|^2) * (F(xn+1) - F(xn) - Gp) */ 332 Gp = (F - F0 - Gp)/PetscSqr(norm); 333 } 334 PetscCall(VecAXPY(G, Gp, Xdiff)); 335 PetscCall(MatMult(S, G , SgF)); /* S*gradF */ 336 337 } else { 338 PetscCall((*dg->Sfunc)(ts, dg->stage_time, x, S, dg->funcCtx)); 339 PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx)); 340 341 PetscCall(MatMult(S, G , SgF)); /* Xdot = S*gradF */ 342 } 343 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 344 dmsave = ts->dm; 345 ts->dm = dm; 346 PetscCall(VecAXPBYPCZ(y, 1, -1, 0, Xdot, SgF)); 347 ts->dm = dmsave; 348 PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 349 350 PetscCall(VecDestroy(&Xp)); 351 PetscCall(VecDestroy(&Xdiff)); 352 PetscCall(VecDestroy(&SgF)); 353 PetscCall(VecDestroy(&G)); 354 PetscCall(MatDestroy(&S)); 355 356 PetscFunctionReturn(0); 357 } 358 359 static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts) 360 { 361 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 362 PetscReal shift = 1/(0.5*ts->time_step); 363 Vec Xdot; 364 DM dm,dmsave; 365 366 PetscFunctionBegin; 367 PetscCall(SNESGetDM(snes, &dm)); 368 /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */ 369 PetscCall(TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot)); 370 371 dmsave = ts->dm; 372 ts->dm = dm; 373 PetscCall(TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE)); 374 ts->dm = dmsave; 375 PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 376 PetscFunctionReturn(0); 377 } 378 379 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) 380 { 381 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 382 383 PetscFunctionBegin; 384 *Sfunc = dg->Sfunc; 385 *Ffunc = dg->Ffunc; 386 *Gfunc = dg->Gfunc; 387 PetscFunctionReturn(0); 388 } 389 390 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) 391 { 392 TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 393 394 PetscFunctionBegin; 395 dg->Sfunc = Sfunc; 396 dg->Ffunc = Ffunc; 397 dg->Gfunc = Gfunc; 398 dg->funcCtx = ctx; 399 PetscFunctionReturn(0); 400 } 401 402 /*MC 403 TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method 404 405 Notes: 406 This is the implicit midpoint rule, with an optional term that guarantees the discrete gradient property. This 407 timestepper applies to systems of the form 408 $ u_t = S(u) grad F(u) 409 where S(u) is a linear operator, and F is a functional of u. 410 411 Level: intermediate 412 413 .seealso: `TSCreate()`, `TSSetType()`, `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()` 414 M*/ 415 PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts) 416 { 417 TS_DiscGrad *th; 418 419 PetscFunctionBegin; 420 PetscCall(PetscCitationsRegister(DGCitation, &DGCite)); 421 ts->ops->reset = TSReset_DiscGrad; 422 ts->ops->destroy = TSDestroy_DiscGrad; 423 ts->ops->view = TSView_DiscGrad; 424 ts->ops->setfromoptions = TSSetFromOptions_DiscGrad; 425 ts->ops->setup = TSSetUp_DiscGrad; 426 ts->ops->step = TSStep_DiscGrad; 427 ts->ops->interpolate = TSInterpolate_DiscGrad; 428 ts->ops->getstages = TSGetStages_DiscGrad; 429 ts->ops->snesfunction = SNESTSFormFunction_DiscGrad; 430 ts->ops->snesjacobian = SNESTSFormJacobian_DiscGrad; 431 ts->default_adapt_type = TSADAPTNONE; 432 433 ts->usessnes = PETSC_TRUE; 434 435 PetscCall(PetscNewLog(ts,&th)); 436 ts->data = (void*)th; 437 438 th->gonzalez = PETSC_FALSE; 439 440 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",TSDiscGradGetFormulation_DiscGrad)); 441 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",TSDiscGradSetFormulation_DiscGrad)); 442 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradIsGonzalez_C",TSDiscGradIsGonzalez_DiscGrad)); 443 PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradUseGonzalez_C",TSDiscGradUseGonzalez_DiscGrad)); 444 PetscFunctionReturn(0); 445 } 446 447 /*@C 448 TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the formulation u_t = S grad F 449 450 Not Collective 451 452 Input Parameter: 453 . ts - timestepping context 454 455 Output Parameters: 456 + Sfunc - constructor for the S matrix from the formulation 457 . Ffunc - functional F from the formulation 458 . Gfunc - constructor for the gradient of F from the formulation 459 - ctx - the user context 460 461 Calling sequence of Sfunc: 462 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 463 464 Calling sequence of Ffunc: 465 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 466 467 Calling sequence of Gfunc: 468 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 469 470 Level: intermediate 471 472 .seealso: `TSDiscGradSetFormulation()` 473 @*/ 474 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) 475 { 476 PetscFunctionBegin; 477 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 478 PetscValidPointer(Sfunc, 2); 479 PetscValidPointer(Ffunc, 3); 480 PetscValidPointer(Gfunc, 4); 481 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)); 482 PetscFunctionReturn(0); 483 } 484 485 /*@C 486 TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u) 487 488 Not Collective 489 490 Input Parameters: 491 + ts - timestepping context 492 . Sfunc - constructor for the S matrix from the formulation 493 . Ffunc - functional F from the formulation 494 - Gfunc - constructor for the gradient of F from the formulation 495 Calling sequence of Sfunc: 496 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 497 498 Calling sequence of Ffunc: 499 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 500 501 Calling sequence of Gfunc: 502 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 503 504 Level: Intermediate 505 506 .seealso: `TSDiscGradGetFormulation()` 507 @*/ 508 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) 509 { 510 PetscFunctionBegin; 511 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 512 PetscValidFunction(Sfunc, 2); 513 PetscValidFunction(Ffunc, 3); 514 PetscValidFunction(Gfunc, 4); 515 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)); 516 PetscFunctionReturn(0); 517 } 518 519 /*@ 520 TSDiscGradIsGonzalez - Checks flag for whether to use additional conservative terms in discrete gradient formulation. 521 522 Not Collective 523 524 Input Parameter: 525 . ts - timestepping context 526 527 Output Parameter: 528 . gonzalez - PETSC_TRUE when using the Gonzalez term 529 530 Level: Advanced 531 532 .seealso: `TSDiscGradUseGonzalez()`, `TSDISCGRAD` 533 @*/ 534 PetscErrorCode TSDiscGradIsGonzalez(TS ts,PetscBool *gonzalez) 535 { 536 PetscFunctionBegin; 537 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 538 PetscValidBoolPointer(gonzalez,2); 539 PetscUseMethod(ts,"TSDiscGradIsGonzalez_C",(TS,PetscBool*),(ts,gonzalez)); 540 PetscFunctionReturn(0); 541 } 542 543 /*@ 544 TSDiscGradUseGonzalez - Sets discrete gradient formulation with or without additional conservative terms. Without flag, the discrete gradients timestepper is just backwards euler 545 546 Not Collective 547 548 Input Parameters: 549 + ts - timestepping context 550 - flg - PETSC_TRUE to use the Gonzalez term 551 552 Options Database: 553 . -ts_discgrad_gonzalez <flg> - use the Gonzalez term for the discrete gradient formulation 554 555 Level: Intermediate 556 557 .seealso: `TSDISCGRAD` 558 @*/ 559 PetscErrorCode TSDiscGradUseGonzalez(TS ts,PetscBool flg) 560 { 561 PetscFunctionBegin; 562 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 563 PetscTryMethod(ts,"TSDiscGradUseGonzalez_C",(TS,PetscBool),(ts,flg)); 564 PetscFunctionReturn(0); 565 } 566