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