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