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