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