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