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 ts->snes_its += its; ts->ksp_its += lits; 188 ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr); 189 ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 190 if (!accept) continue; 191 ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr); 192 /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 193 ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr); 194 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 195 ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr); 196 ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 197 198 if (accept) { 199 /* ignore next_scheme for now */ 200 ts->ptime += ts->time_step; 201 ts->time_step = next_time_step; 202 ts->steps++; 203 th->status = TS_STEP_COMPLETE; 204 } else { /* Roll back the current step */ 205 ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 206 ts->time_step = next_time_step; 207 th->status = TS_STEP_INCOMPLETE; 208 } 209 } 210 PetscFunctionReturn(0); 211 } 212 213 #undef __FUNCT__ 214 #define __FUNCT__ "TSInterpolate_Theta" 215 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 216 { 217 TS_Theta *th = (TS_Theta*)ts->data; 218 PetscReal alpha = t - ts->ptime; 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 223 if (th->endpoint) alpha *= th->Theta; 224 ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 225 PetscFunctionReturn(0); 226 } 227 228 /*------------------------------------------------------------*/ 229 #undef __FUNCT__ 230 #define __FUNCT__ "TSReset_Theta" 231 static PetscErrorCode TSReset_Theta(TS ts) 232 { 233 TS_Theta *th = (TS_Theta*)ts->data; 234 PetscErrorCode ierr; 235 236 PetscFunctionBegin; 237 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 238 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 239 ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 240 ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 241 PetscFunctionReturn(0); 242 } 243 244 #undef __FUNCT__ 245 #define __FUNCT__ "TSDestroy_Theta" 246 static PetscErrorCode TSDestroy_Theta(TS ts) 247 { 248 PetscErrorCode ierr; 249 250 PetscFunctionBegin; 251 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 252 ierr = PetscFree(ts->data);CHKERRQ(ierr); 253 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",NULL);CHKERRQ(ierr); 254 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",NULL);CHKERRQ(ierr); 255 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",NULL);CHKERRQ(ierr); 256 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",NULL);CHKERRQ(ierr); 257 PetscFunctionReturn(0); 258 } 259 260 /* 261 This defines the nonlinear equation that is to be solved with SNES 262 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 263 */ 264 #undef __FUNCT__ 265 #define __FUNCT__ "SNESTSFormFunction_Theta" 266 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 267 { 268 TS_Theta *th = (TS_Theta*)ts->data; 269 PetscErrorCode ierr; 270 Vec X0,Xdot; 271 DM dm,dmsave; 272 PetscReal shift = 1./(th->Theta*ts->time_step); 273 274 PetscFunctionBegin; 275 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 276 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 277 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 278 ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 279 280 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 281 dmsave = ts->dm; 282 ts->dm = dm; 283 ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 284 ts->dm = dmsave; 285 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "SNESTSFormJacobian_Theta" 291 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 292 { 293 TS_Theta *th = (TS_Theta*)ts->data; 294 PetscErrorCode ierr; 295 Vec Xdot; 296 DM dm,dmsave; 297 PetscReal shift = 1./(th->Theta*ts->time_step); 298 299 PetscFunctionBegin; 300 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 301 302 /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 303 ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 304 305 dmsave = ts->dm; 306 ts->dm = dm; 307 ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 308 ts->dm = dmsave; 309 ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 310 PetscFunctionReturn(0); 311 } 312 313 #undef __FUNCT__ 314 #define __FUNCT__ "TSSetUp_Theta" 315 static PetscErrorCode TSSetUp_Theta(TS ts) 316 { 317 TS_Theta *th = (TS_Theta*)ts->data; 318 PetscErrorCode ierr; 319 SNES snes; 320 DM dm; 321 322 PetscFunctionBegin; 323 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 324 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 325 ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 326 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 327 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 328 if (dm) { 329 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 330 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 331 } 332 if (th->Theta == 0.5 && th->endpoint) th->order = 2; 333 else th->order = 1; 334 335 if (!th->adapt) { 336 TSAdapt adapt; 337 ierr = TSAdaptDestroy(&ts->adapt);CHKERRQ(ierr); 338 ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr); 339 ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr); 340 } 341 PetscFunctionReturn(0); 342 } 343 /*------------------------------------------------------------*/ 344 345 #undef __FUNCT__ 346 #define __FUNCT__ "TSSetFromOptions_Theta" 347 static PetscErrorCode TSSetFromOptions_Theta(TS ts) 348 { 349 TS_Theta *th = (TS_Theta*)ts->data; 350 PetscErrorCode ierr; 351 352 PetscFunctionBegin; 353 ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 354 { 355 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 356 ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 357 ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr); 358 ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr); 359 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 360 } 361 ierr = PetscOptionsTail();CHKERRQ(ierr); 362 PetscFunctionReturn(0); 363 } 364 365 #undef __FUNCT__ 366 #define __FUNCT__ "TSView_Theta" 367 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 368 { 369 TS_Theta *th = (TS_Theta*)ts->data; 370 PetscBool iascii; 371 PetscErrorCode ierr; 372 373 PetscFunctionBegin; 374 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 375 if (iascii) { 376 ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 377 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 378 } 379 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 380 PetscFunctionReturn(0); 381 } 382 383 EXTERN_C_BEGIN 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 EXTERN_C_END 429 430 #if defined(PETSC_HAVE_COMPLEX) 431 #undef __FUNCT__ 432 #define __FUNCT__ "TSComputeLinearStability_Theta" 433 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 434 { 435 PetscComplex z = xr + xi*PETSC_i,f; 436 TS_Theta *th = (TS_Theta*)ts->data; 437 const PetscReal one = 1.0; 438 439 PetscFunctionBegin; 440 f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 441 *yr = PetscRealPartComplex(f); 442 *yi = PetscImaginaryPartComplex(f); 443 PetscFunctionReturn(0); 444 } 445 #endif 446 447 448 /* ------------------------------------------------------------ */ 449 /*MC 450 TSTHETA - DAE solver using the implicit Theta method 451 452 Level: beginner 453 454 Options Database: 455 -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 456 -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 457 -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 458 459 Notes: 460 $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 461 $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 462 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 463 464 465 466 This method can be applied to DAE. 467 468 This method is cast as a 1-stage implicit Runge-Kutta method. 469 470 .vb 471 Theta | Theta 472 ------------- 473 | 1 474 .ve 475 476 For the default Theta=0.5, this is also known as the implicit midpoint rule. 477 478 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 479 480 .vb 481 0 | 0 0 482 1 | 1-Theta Theta 483 ------------------- 484 | 1-Theta Theta 485 .ve 486 487 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 488 489 To apply a diagonally implicit RK method to DAE, the stage formula 490 491 $ Y_i = X + h sum_j a_ij Y'_j 492 493 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 494 495 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 496 497 M*/ 498 EXTERN_C_BEGIN 499 #undef __FUNCT__ 500 #define __FUNCT__ "TSCreate_Theta" 501 PetscErrorCode TSCreate_Theta(TS ts) 502 { 503 TS_Theta *th; 504 PetscErrorCode ierr; 505 506 PetscFunctionBegin; 507 ts->ops->reset = TSReset_Theta; 508 ts->ops->destroy = TSDestroy_Theta; 509 ts->ops->view = TSView_Theta; 510 ts->ops->setup = TSSetUp_Theta; 511 ts->ops->step = TSStep_Theta; 512 ts->ops->interpolate = TSInterpolate_Theta; 513 ts->ops->evaluatestep = TSEvaluateStep_Theta; 514 ts->ops->setfromoptions = TSSetFromOptions_Theta; 515 ts->ops->snesfunction = SNESTSFormFunction_Theta; 516 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 517 #if defined(PETSC_HAVE_COMPLEX) 518 ts->ops->linearstability = TSComputeLinearStability_Theta; 519 #endif 520 521 ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 522 ts->data = (void*)th; 523 524 th->extrapolate = PETSC_FALSE; 525 th->Theta = 0.5; 526 th->ccfl = 1.0; 527 th->adapt = PETSC_FALSE; 528 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 529 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 530 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 531 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 532 PetscFunctionReturn(0); 533 } 534 EXTERN_C_END 535 536 #undef __FUNCT__ 537 #define __FUNCT__ "TSThetaGetTheta" 538 /*@ 539 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 540 541 Not Collective 542 543 Input Parameter: 544 . ts - timestepping context 545 546 Output Parameter: 547 . theta - stage abscissa 548 549 Note: 550 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 551 552 Level: Advanced 553 554 .seealso: TSThetaSetTheta() 555 @*/ 556 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 557 { 558 PetscErrorCode ierr; 559 560 PetscFunctionBegin; 561 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 562 PetscValidPointer(theta,2); 563 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 564 PetscFunctionReturn(0); 565 } 566 567 #undef __FUNCT__ 568 #define __FUNCT__ "TSThetaSetTheta" 569 /*@ 570 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 571 572 Not Collective 573 574 Input Parameter: 575 + ts - timestepping context 576 - theta - stage abscissa 577 578 Options Database: 579 . -ts_theta_theta <theta> 580 581 Level: Intermediate 582 583 .seealso: TSThetaGetTheta() 584 @*/ 585 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 586 { 587 PetscErrorCode ierr; 588 589 PetscFunctionBegin; 590 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 591 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 592 PetscFunctionReturn(0); 593 } 594 595 #undef __FUNCT__ 596 #define __FUNCT__ "TSThetaGetEndpoint" 597 /*@ 598 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 599 600 Not Collective 601 602 Input Parameter: 603 . ts - timestepping context 604 605 Output Parameter: 606 . endpoint - PETSC_TRUE when using the endpoint variant 607 608 Level: Advanced 609 610 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 611 @*/ 612 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 613 { 614 PetscErrorCode ierr; 615 616 PetscFunctionBegin; 617 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 618 PetscValidPointer(endpoint,2); 619 ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 620 PetscFunctionReturn(0); 621 } 622 623 #undef __FUNCT__ 624 #define __FUNCT__ "TSThetaSetEndpoint" 625 /*@ 626 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 627 628 Not Collective 629 630 Input Parameter: 631 + ts - timestepping context 632 - flg - PETSC_TRUE to use the endpoint variant 633 634 Options Database: 635 . -ts_theta_endpoint <flg> 636 637 Level: Intermediate 638 639 .seealso: TSTHETA, TSCN 640 @*/ 641 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 642 { 643 PetscErrorCode ierr; 644 645 PetscFunctionBegin; 646 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 647 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 648 PetscFunctionReturn(0); 649 } 650 651 /* 652 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 653 * The creation functions for these specializations are below. 654 */ 655 656 #undef __FUNCT__ 657 #define __FUNCT__ "TSView_BEuler" 658 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 659 { 660 PetscErrorCode ierr; 661 662 PetscFunctionBegin; 663 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 664 PetscFunctionReturn(0); 665 } 666 667 /*MC 668 TSBEULER - ODE solver using the implicit backward Euler method 669 670 Level: beginner 671 672 Notes: 673 TSBEULER is equivalent to TSTHETA with Theta=1.0 674 675 $ -ts_type theta -ts_theta_theta 1. 676 677 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 678 679 M*/ 680 EXTERN_C_BEGIN 681 #undef __FUNCT__ 682 #define __FUNCT__ "TSCreate_BEuler" 683 PetscErrorCode TSCreate_BEuler(TS ts) 684 { 685 PetscErrorCode ierr; 686 687 PetscFunctionBegin; 688 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 689 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 690 ts->ops->view = TSView_BEuler; 691 PetscFunctionReturn(0); 692 } 693 EXTERN_C_END 694 695 #undef __FUNCT__ 696 #define __FUNCT__ "TSView_CN" 697 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 698 { 699 PetscErrorCode ierr; 700 701 PetscFunctionBegin; 702 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 703 PetscFunctionReturn(0); 704 } 705 706 /*MC 707 TSCN - ODE solver using the implicit Crank-Nicolson method. 708 709 Level: beginner 710 711 Notes: 712 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 713 714 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 715 716 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 717 718 M*/ 719 EXTERN_C_BEGIN 720 #undef __FUNCT__ 721 #define __FUNCT__ "TSCreate_CN" 722 PetscErrorCode TSCreate_CN(TS ts) 723 { 724 PetscErrorCode ierr; 725 726 PetscFunctionBegin; 727 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 728 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 729 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 730 ts->ops->view = TSView_CN; 731 PetscFunctionReturn(0); 732 } 733 EXTERN_C_END 734