1 /* 2 Code for timestepping with implicit Theta method 3 */ 4 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 5 #include <petscsnesfas.h> 6 7 typedef struct { 8 Vec X,Xdot; /* Storage for one stage */ 9 Vec affine; /* Affine vector needed for residual at beginning of step */ 10 PetscBool extrapolate; 11 PetscBool endpoint; 12 PetscReal Theta; 13 PetscReal shift; 14 PetscReal stage_time; 15 } TS_Theta; 16 17 #undef __FUNCT__ 18 #define __FUNCT__ "TSThetaGetX0AndXdot" 19 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 20 { 21 TS_Theta *th = (TS_Theta*)ts->data; 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 if (X0) { 26 if (dm && dm != ts->dm) { 27 ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 28 } else *X0 = ts->vec_sol; 29 } 30 if (Xdot) { 31 if (dm && dm != ts->dm) { 32 ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 33 } else *Xdot = th->Xdot; 34 } 35 PetscFunctionReturn(0); 36 } 37 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "TSThetaRestoreX0AndXdot" 41 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 42 { 43 PetscErrorCode ierr; 44 45 PetscFunctionBegin; 46 if (X0) { 47 if (dm && dm != ts->dm) { 48 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 49 } 50 } 51 if (Xdot) { 52 if (dm && dm != ts->dm) { 53 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 54 } 55 } 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "DMCoarsenHook_TSTheta" 61 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 62 { 63 64 PetscFunctionBegin; 65 PetscFunctionReturn(0); 66 } 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "DMRestrictHook_TSTheta" 70 static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 71 { 72 TS ts = (TS)ctx; 73 PetscErrorCode ierr; 74 Vec X0,Xdot,X0_c,Xdot_c; 75 76 PetscFunctionBegin; 77 ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 78 ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 79 ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 80 ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 81 ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 82 ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 83 ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 84 ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 85 PetscFunctionReturn(0); 86 } 87 88 #undef __FUNCT__ 89 #define __FUNCT__ "TSStep_Theta" 90 static PetscErrorCode TSStep_Theta(TS ts) 91 { 92 TS_Theta *th = (TS_Theta*)ts->data; 93 PetscInt its,lits,reject; 94 PetscReal next_time_step; 95 SNESConvergedReason snesreason; 96 PetscErrorCode ierr; 97 98 PetscFunctionBegin; 99 if (ts->time_steps_since_decrease > 3 && ts->time_step < ts->time_step_orig) { 100 /* smaller time step has worked successfully for three time-steps, try increasing time step*/ 101 ts->time_step = 2.0*ts->time_step; 102 ts->time_steps_since_decrease = 0; /* don't want to increase time step two time steps in a row */ 103 } 104 for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 105 next_time_step = ts->time_step; 106 th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 107 th->shift = 1./(th->Theta*ts->time_step); 108 ierr = TSPreStep(ts);CHKERRQ(ierr); 109 ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 110 111 if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 112 ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 113 if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 114 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 115 ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 116 } 117 if (th->extrapolate) { 118 ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 119 } else { 120 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 121 } 122 ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 123 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 124 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 125 ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 126 ts->snes_its += its; ts->ksp_its += lits; 127 if (its < 10) ts->time_steps_since_decrease++; 128 else ts->time_steps_since_decrease = 0; 129 if (snesreason > 0) break; 130 ierr = PetscInfo3(ts,"Step=%D, Cutting time-step from %g to %g\n",ts->steps,(double)ts->time_step,(double).5*ts->time_step);CHKERRQ(ierr); 131 ts->time_step = .5*ts->time_step; 132 ts->time_steps_since_decrease = 0; 133 } 134 if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 135 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 136 ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 137 PetscFunctionReturn(0); 138 } 139 if (th->endpoint) { 140 ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 141 } else { 142 ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 143 ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 144 } 145 ts->ptime += ts->time_step; 146 ts->time_step = next_time_step; 147 ts->steps++; 148 PetscFunctionReturn(0); 149 } 150 151 #undef __FUNCT__ 152 #define __FUNCT__ "TSInterpolate_Theta" 153 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 154 { 155 TS_Theta *th = (TS_Theta*)ts->data; 156 PetscReal alpha = t - ts->ptime; 157 PetscErrorCode ierr; 158 159 PetscFunctionBegin; 160 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 161 if (th->endpoint) alpha *= th->Theta; 162 ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 163 PetscFunctionReturn(0); 164 } 165 166 /*------------------------------------------------------------*/ 167 #undef __FUNCT__ 168 #define __FUNCT__ "TSReset_Theta" 169 static PetscErrorCode TSReset_Theta(TS ts) 170 { 171 TS_Theta *th = (TS_Theta*)ts->data; 172 PetscErrorCode ierr; 173 174 PetscFunctionBegin; 175 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 176 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 177 ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "TSDestroy_Theta" 183 static PetscErrorCode TSDestroy_Theta(TS ts) 184 { 185 PetscErrorCode ierr; 186 187 PetscFunctionBegin; 188 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 189 ierr = PetscFree(ts->data);CHKERRQ(ierr); 190 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 191 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 192 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 193 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 /* 198 This defines the nonlinear equation that is to be solved with SNES 199 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 200 */ 201 #undef __FUNCT__ 202 #define __FUNCT__ "SNESTSFormFunction_Theta" 203 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 204 { 205 TS_Theta *th = (TS_Theta*)ts->data; 206 PetscErrorCode ierr; 207 Vec X0,Xdot; 208 DM dm,dmsave; 209 210 PetscFunctionBegin; 211 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 212 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 213 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 214 ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr); 215 216 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 217 dmsave = ts->dm; 218 ts->dm = dm; 219 ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 220 ts->dm = dmsave; 221 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 222 PetscFunctionReturn(0); 223 } 224 225 #undef __FUNCT__ 226 #define __FUNCT__ "SNESTSFormJacobian_Theta" 227 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 228 { 229 TS_Theta *th = (TS_Theta*)ts->data; 230 PetscErrorCode ierr; 231 Vec Xdot; 232 DM dm,dmsave; 233 234 PetscFunctionBegin; 235 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 236 237 /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 238 ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 239 240 dmsave = ts->dm; 241 ts->dm = dm; 242 ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 243 ts->dm = dmsave; 244 ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "TSSetUp_Theta" 250 static PetscErrorCode TSSetUp_Theta(TS ts) 251 { 252 TS_Theta *th = (TS_Theta*)ts->data; 253 PetscErrorCode ierr; 254 SNES snes; 255 DM dm; 256 257 PetscFunctionBegin; 258 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 259 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 260 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 261 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 262 if (dm) { 263 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 264 } 265 PetscFunctionReturn(0); 266 } 267 /*------------------------------------------------------------*/ 268 269 #undef __FUNCT__ 270 #define __FUNCT__ "TSSetFromOptions_Theta" 271 static PetscErrorCode TSSetFromOptions_Theta(TS ts) 272 { 273 TS_Theta *th = (TS_Theta*)ts->data; 274 PetscErrorCode ierr; 275 276 PetscFunctionBegin; 277 ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 278 { 279 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 280 ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 281 ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);CHKERRQ(ierr); 282 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 283 } 284 ierr = PetscOptionsTail();CHKERRQ(ierr); 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "TSView_Theta" 290 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 291 { 292 TS_Theta *th = (TS_Theta*)ts->data; 293 PetscBool iascii; 294 PetscErrorCode ierr; 295 296 PetscFunctionBegin; 297 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 298 if (iascii) { 299 ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 300 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 301 } 302 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 EXTERN_C_BEGIN 307 #undef __FUNCT__ 308 #define __FUNCT__ "TSThetaGetTheta_Theta" 309 PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 310 { 311 TS_Theta *th = (TS_Theta*)ts->data; 312 313 PetscFunctionBegin; 314 *theta = th->Theta; 315 PetscFunctionReturn(0); 316 } 317 318 #undef __FUNCT__ 319 #define __FUNCT__ "TSThetaSetTheta_Theta" 320 PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 321 { 322 TS_Theta *th = (TS_Theta*)ts->data; 323 324 PetscFunctionBegin; 325 if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 326 th->Theta = theta; 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "TSThetaGetEndpoint_Theta" 332 PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 333 { 334 TS_Theta *th = (TS_Theta*)ts->data; 335 336 PetscFunctionBegin; 337 *endpoint = th->endpoint; 338 PetscFunctionReturn(0); 339 } 340 341 #undef __FUNCT__ 342 #define __FUNCT__ "TSThetaSetEndpoint_Theta" 343 PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 344 { 345 TS_Theta *th = (TS_Theta*)ts->data; 346 347 PetscFunctionBegin; 348 th->endpoint = flg; 349 PetscFunctionReturn(0); 350 } 351 EXTERN_C_END 352 353 /* ------------------------------------------------------------ */ 354 /*MC 355 TSTHETA - DAE solver using the implicit Theta method 356 357 Level: beginner 358 359 Options Database: 360 -ts_theta_theta <Theta> - Location of stage (0<Theta<=1); Theta = 1.0 ( 361 -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 362 -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method 363 364 Notes: 365 $ -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER) 366 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 367 368 369 370 This method can be applied to DAE. 371 372 This method is cast as a 1-stage implicit Runge-Kutta method. 373 374 .vb 375 Theta | Theta 376 ------------- 377 | 1 378 .ve 379 380 For the default Theta=0.5, this is also known as the implicit midpoint rule. 381 382 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 383 384 .vb 385 0 | 0 0 386 1 | 1-Theta Theta 387 ------------------- 388 | 1-Theta Theta 389 .ve 390 391 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 392 393 To apply a diagonally implicit RK method to DAE, the stage formula 394 395 $ Y_i = X + h sum_j a_ij Y'_j 396 397 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 398 399 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 400 401 M*/ 402 EXTERN_C_BEGIN 403 #undef __FUNCT__ 404 #define __FUNCT__ "TSCreate_Theta" 405 PetscErrorCode TSCreate_Theta(TS ts) 406 { 407 TS_Theta *th; 408 PetscErrorCode ierr; 409 410 PetscFunctionBegin; 411 ts->ops->reset = TSReset_Theta; 412 ts->ops->destroy = TSDestroy_Theta; 413 ts->ops->view = TSView_Theta; 414 ts->ops->setup = TSSetUp_Theta; 415 ts->ops->step = TSStep_Theta; 416 ts->ops->interpolate = TSInterpolate_Theta; 417 ts->ops->setfromoptions = TSSetFromOptions_Theta; 418 ts->ops->snesfunction = SNESTSFormFunction_Theta; 419 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 420 421 ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 422 ts->data = (void*)th; 423 424 th->extrapolate = PETSC_FALSE; 425 th->Theta = 0.5; 426 427 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 428 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 429 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 430 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 431 PetscFunctionReturn(0); 432 } 433 EXTERN_C_END 434 435 #undef __FUNCT__ 436 #define __FUNCT__ "TSThetaGetTheta" 437 /*@ 438 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 439 440 Not Collective 441 442 Input Parameter: 443 . ts - timestepping context 444 445 Output Parameter: 446 . theta - stage abscissa 447 448 Note: 449 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 450 451 Level: Advanced 452 453 .seealso: TSThetaSetTheta() 454 @*/ 455 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 456 { 457 PetscErrorCode ierr; 458 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 461 PetscValidPointer(theta,2); 462 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNCT__ 467 #define __FUNCT__ "TSThetaSetTheta" 468 /*@ 469 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 470 471 Not Collective 472 473 Input Parameter: 474 + ts - timestepping context 475 - theta - stage abscissa 476 477 Options Database: 478 . -ts_theta_theta <theta> 479 480 Level: Intermediate 481 482 .seealso: TSThetaGetTheta() 483 @*/ 484 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 485 { 486 PetscErrorCode ierr; 487 488 PetscFunctionBegin; 489 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 490 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "TSThetaGetEndpoint" 496 /*@ 497 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 498 499 Not Collective 500 501 Input Parameter: 502 . ts - timestepping context 503 504 Output Parameter: 505 . endpoint - PETSC_TRUE when using the endpoint variant 506 507 Level: Advanced 508 509 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 510 @*/ 511 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 512 { 513 PetscErrorCode ierr; 514 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 517 PetscValidPointer(endpoint,2); 518 ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 519 PetscFunctionReturn(0); 520 } 521 522 #undef __FUNCT__ 523 #define __FUNCT__ "TSThetaSetEndpoint" 524 /*@ 525 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 526 527 Not Collective 528 529 Input Parameter: 530 + ts - timestepping context 531 - flg - PETSC_TRUE to use the endpoint variant 532 533 Options Database: 534 . -ts_theta_endpoint <flg> 535 536 Level: Intermediate 537 538 .seealso: TSTHETA, TSCN 539 @*/ 540 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 541 { 542 PetscErrorCode ierr; 543 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 546 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 547 PetscFunctionReturn(0); 548 } 549 550 /* 551 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 552 * The creation functions for these specializations are below. 553 */ 554 555 #undef __FUNCT__ 556 #define __FUNCT__ "TSView_BEuler" 557 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 558 { 559 PetscErrorCode ierr; 560 561 PetscFunctionBegin; 562 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 563 PetscFunctionReturn(0); 564 } 565 566 /*MC 567 TSBEULER - ODE solver using the implicit backward Euler method 568 569 Level: beginner 570 571 Notes: 572 TSCN is equivalent to TSTHETA with Theta=1.0 573 574 $ -ts_type theta -ts_theta_theta 1. 575 576 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 577 578 M*/ 579 EXTERN_C_BEGIN 580 #undef __FUNCT__ 581 #define __FUNCT__ "TSCreate_BEuler" 582 PetscErrorCode TSCreate_BEuler(TS ts) 583 { 584 PetscErrorCode ierr; 585 586 PetscFunctionBegin; 587 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 588 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 589 ts->ops->view = TSView_BEuler; 590 PetscFunctionReturn(0); 591 } 592 EXTERN_C_END 593 594 #undef __FUNCT__ 595 #define __FUNCT__ "TSView_CN" 596 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 597 { 598 PetscErrorCode ierr; 599 600 PetscFunctionBegin; 601 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 602 PetscFunctionReturn(0); 603 } 604 605 /*MC 606 TSCN - ODE solver using the implicit Crank-Nicolson method. 607 608 Level: beginner 609 610 Notes: 611 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 612 613 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 614 615 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 616 617 M*/ 618 EXTERN_C_BEGIN 619 #undef __FUNCT__ 620 #define __FUNCT__ "TSCreate_CN" 621 PetscErrorCode TSCreate_CN(TS ts) 622 { 623 PetscErrorCode ierr; 624 625 PetscFunctionBegin; 626 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 627 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 628 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 629 ts->ops->view = TSView_CN; 630 PetscFunctionReturn(0); 631 } 632 EXTERN_C_END 633