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 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 = 0.0; 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 PetscReal shift = 1./(th->Theta*ts->time_step); 106 next_time_step = ts->time_step; 107 th->stage_time = ts->ptime + (th->endpoint ? 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./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 PetscReal shift = 1./(th->Theta*ts->time_step); 143 ierr = VecAXPBYPCZ(th->Xdot,-shift,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 PetscReal shift = 1./(th->Theta*ts->time_step); 211 212 PetscFunctionBegin; 213 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 214 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 215 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 216 ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 217 218 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 219 dmsave = ts->dm; 220 ts->dm = dm; 221 ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 222 ts->dm = dmsave; 223 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 224 PetscFunctionReturn(0); 225 } 226 227 #undef __FUNCT__ 228 #define __FUNCT__ "SNESTSFormJacobian_Theta" 229 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 230 { 231 TS_Theta *th = (TS_Theta*)ts->data; 232 PetscErrorCode ierr; 233 Vec Xdot; 234 DM dm,dmsave; 235 PetscReal shift = 1./(th->Theta*ts->time_step); 236 237 PetscFunctionBegin; 238 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 239 240 /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 241 ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 242 243 dmsave = ts->dm; 244 ts->dm = dm; 245 ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 246 ts->dm = dmsave; 247 ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "TSSetUp_Theta" 253 static PetscErrorCode TSSetUp_Theta(TS ts) 254 { 255 TS_Theta *th = (TS_Theta*)ts->data; 256 PetscErrorCode ierr; 257 SNES snes; 258 DM dm; 259 260 PetscFunctionBegin; 261 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 262 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 263 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 264 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 265 if (dm) { 266 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 267 } 268 PetscFunctionReturn(0); 269 } 270 /*------------------------------------------------------------*/ 271 272 #undef __FUNCT__ 273 #define __FUNCT__ "TSSetFromOptions_Theta" 274 static PetscErrorCode TSSetFromOptions_Theta(TS ts) 275 { 276 TS_Theta *th = (TS_Theta*)ts->data; 277 PetscErrorCode ierr; 278 279 PetscFunctionBegin; 280 ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 281 { 282 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 283 ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 284 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); 285 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 286 } 287 ierr = PetscOptionsTail();CHKERRQ(ierr); 288 PetscFunctionReturn(0); 289 } 290 291 #undef __FUNCT__ 292 #define __FUNCT__ "TSView_Theta" 293 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 294 { 295 TS_Theta *th = (TS_Theta*)ts->data; 296 PetscBool iascii; 297 PetscErrorCode ierr; 298 299 PetscFunctionBegin; 300 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 301 if (iascii) { 302 ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 303 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 304 } 305 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 306 PetscFunctionReturn(0); 307 } 308 309 EXTERN_C_BEGIN 310 #undef __FUNCT__ 311 #define __FUNCT__ "TSThetaGetTheta_Theta" 312 PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 313 { 314 TS_Theta *th = (TS_Theta*)ts->data; 315 316 PetscFunctionBegin; 317 *theta = th->Theta; 318 PetscFunctionReturn(0); 319 } 320 321 #undef __FUNCT__ 322 #define __FUNCT__ "TSThetaSetTheta_Theta" 323 PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 324 { 325 TS_Theta *th = (TS_Theta*)ts->data; 326 327 PetscFunctionBegin; 328 if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 329 th->Theta = theta; 330 PetscFunctionReturn(0); 331 } 332 333 #undef __FUNCT__ 334 #define __FUNCT__ "TSThetaGetEndpoint_Theta" 335 PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 336 { 337 TS_Theta *th = (TS_Theta*)ts->data; 338 339 PetscFunctionBegin; 340 *endpoint = th->endpoint; 341 PetscFunctionReturn(0); 342 } 343 344 #undef __FUNCT__ 345 #define __FUNCT__ "TSThetaSetEndpoint_Theta" 346 PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 347 { 348 TS_Theta *th = (TS_Theta*)ts->data; 349 350 PetscFunctionBegin; 351 th->endpoint = flg; 352 PetscFunctionReturn(0); 353 } 354 EXTERN_C_END 355 356 #if defined(PETSC_HAVE_COMPLEX) 357 #undef __FUNCT__ 358 #define __FUNCT__ "TSComputeLinearStability_Theta" 359 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 360 { 361 PetscComplex z = xr + xi*PETSC_i,f; 362 TS_Theta *th = (TS_Theta*)ts->data; 363 364 PetscFunctionBegin; 365 f = (1.0 + (1.0 - th->Theta)*z)/(1.0 - th->Theta*z); 366 *yr = PetscRealPartComplex(f); 367 *yi = PetscImaginaryPartComplex(f); 368 PetscFunctionReturn(0); 369 } 370 #endif 371 372 373 /* ------------------------------------------------------------ */ 374 /*MC 375 TSTHETA - DAE solver using the implicit Theta method 376 377 Level: beginner 378 379 Options Database: 380 -ts_theta_theta <Theta> - Location of stage (0<Theta<=1); Theta = 1.0 ( 381 -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 382 -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method 383 384 Notes: 385 $ -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER) 386 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 387 388 389 390 This method can be applied to DAE. 391 392 This method is cast as a 1-stage implicit Runge-Kutta method. 393 394 .vb 395 Theta | Theta 396 ------------- 397 | 1 398 .ve 399 400 For the default Theta=0.5, this is also known as the implicit midpoint rule. 401 402 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 403 404 .vb 405 0 | 0 0 406 1 | 1-Theta Theta 407 ------------------- 408 | 1-Theta Theta 409 .ve 410 411 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 412 413 To apply a diagonally implicit RK method to DAE, the stage formula 414 415 $ Y_i = X + h sum_j a_ij Y'_j 416 417 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 418 419 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 420 421 M*/ 422 EXTERN_C_BEGIN 423 #undef __FUNCT__ 424 #define __FUNCT__ "TSCreate_Theta" 425 PetscErrorCode TSCreate_Theta(TS ts) 426 { 427 TS_Theta *th; 428 PetscErrorCode ierr; 429 430 PetscFunctionBegin; 431 ts->ops->reset = TSReset_Theta; 432 ts->ops->destroy = TSDestroy_Theta; 433 ts->ops->view = TSView_Theta; 434 ts->ops->setup = TSSetUp_Theta; 435 ts->ops->step = TSStep_Theta; 436 ts->ops->interpolate = TSInterpolate_Theta; 437 ts->ops->setfromoptions = TSSetFromOptions_Theta; 438 ts->ops->snesfunction = SNESTSFormFunction_Theta; 439 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 440 #if defined(PETSC_HAVE_COMPLEX) 441 ts->ops->linearstability = TSComputeLinearStability_Theta; 442 #endif 443 444 ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 445 ts->data = (void*)th; 446 447 th->extrapolate = PETSC_FALSE; 448 th->Theta = 0.5; 449 450 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 451 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 452 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 453 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 454 PetscFunctionReturn(0); 455 } 456 EXTERN_C_END 457 458 #undef __FUNCT__ 459 #define __FUNCT__ "TSThetaGetTheta" 460 /*@ 461 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 462 463 Not Collective 464 465 Input Parameter: 466 . ts - timestepping context 467 468 Output Parameter: 469 . theta - stage abscissa 470 471 Note: 472 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 473 474 Level: Advanced 475 476 .seealso: TSThetaSetTheta() 477 @*/ 478 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 479 { 480 PetscErrorCode ierr; 481 482 PetscFunctionBegin; 483 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 484 PetscValidPointer(theta,2); 485 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 486 PetscFunctionReturn(0); 487 } 488 489 #undef __FUNCT__ 490 #define __FUNCT__ "TSThetaSetTheta" 491 /*@ 492 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 493 494 Not Collective 495 496 Input Parameter: 497 + ts - timestepping context 498 - theta - stage abscissa 499 500 Options Database: 501 . -ts_theta_theta <theta> 502 503 Level: Intermediate 504 505 .seealso: TSThetaGetTheta() 506 @*/ 507 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 508 { 509 PetscErrorCode ierr; 510 511 PetscFunctionBegin; 512 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 513 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 514 PetscFunctionReturn(0); 515 } 516 517 #undef __FUNCT__ 518 #define __FUNCT__ "TSThetaGetEndpoint" 519 /*@ 520 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 521 522 Not Collective 523 524 Input Parameter: 525 . ts - timestepping context 526 527 Output Parameter: 528 . endpoint - PETSC_TRUE when using the endpoint variant 529 530 Level: Advanced 531 532 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 533 @*/ 534 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 535 { 536 PetscErrorCode ierr; 537 538 PetscFunctionBegin; 539 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 540 PetscValidPointer(endpoint,2); 541 ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 542 PetscFunctionReturn(0); 543 } 544 545 #undef __FUNCT__ 546 #define __FUNCT__ "TSThetaSetEndpoint" 547 /*@ 548 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 549 550 Not Collective 551 552 Input Parameter: 553 + ts - timestepping context 554 - flg - PETSC_TRUE to use the endpoint variant 555 556 Options Database: 557 . -ts_theta_endpoint <flg> 558 559 Level: Intermediate 560 561 .seealso: TSTHETA, TSCN 562 @*/ 563 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 564 { 565 PetscErrorCode ierr; 566 567 PetscFunctionBegin; 568 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 569 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 570 PetscFunctionReturn(0); 571 } 572 573 /* 574 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 575 * The creation functions for these specializations are below. 576 */ 577 578 #undef __FUNCT__ 579 #define __FUNCT__ "TSView_BEuler" 580 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 581 { 582 PetscErrorCode ierr; 583 584 PetscFunctionBegin; 585 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 586 PetscFunctionReturn(0); 587 } 588 589 /*MC 590 TSBEULER - ODE solver using the implicit backward Euler method 591 592 Level: beginner 593 594 Notes: 595 TSCN is equivalent to TSTHETA with Theta=1.0 596 597 $ -ts_type theta -ts_theta_theta 1. 598 599 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 600 601 M*/ 602 EXTERN_C_BEGIN 603 #undef __FUNCT__ 604 #define __FUNCT__ "TSCreate_BEuler" 605 PetscErrorCode TSCreate_BEuler(TS ts) 606 { 607 PetscErrorCode ierr; 608 609 PetscFunctionBegin; 610 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 611 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 612 ts->ops->view = TSView_BEuler; 613 PetscFunctionReturn(0); 614 } 615 EXTERN_C_END 616 617 #undef __FUNCT__ 618 #define __FUNCT__ "TSView_CN" 619 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 620 { 621 PetscErrorCode ierr; 622 623 PetscFunctionBegin; 624 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 625 PetscFunctionReturn(0); 626 } 627 628 /*MC 629 TSCN - ODE solver using the implicit Crank-Nicolson method. 630 631 Level: beginner 632 633 Notes: 634 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 635 636 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 637 638 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 639 640 M*/ 641 EXTERN_C_BEGIN 642 #undef __FUNCT__ 643 #define __FUNCT__ "TSCreate_CN" 644 PetscErrorCode TSCreate_CN(TS ts) 645 { 646 PetscErrorCode ierr; 647 648 PetscFunctionBegin; 649 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 650 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 651 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 652 ts->ops->view = TSView_CN; 653 PetscFunctionReturn(0); 654 } 655 EXTERN_C_END 656