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.) { /* two-stage case */ 280 shift = -1./((th->Theta-1.)*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 287 if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 288 ierr = TSRHSJacobianP(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 289 for (nadj=0; nadj<ts->numberadjs; nadj++) { 290 ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr); 291 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecDeltaMu[nadj]);CHKERRQ(ierr); 292 } 293 ierr = TSRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 294 for (nadj=0; nadj<ts->numberadjs; nadj++) { 295 ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr); 296 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecDeltaMu[nadj]);CHKERRQ(ierr); 297 } 298 } 299 }else { /* one-stage case */ 300 shift = 0.0; 301 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 302 /* When th->endpoint is true and th->Theta==1 (beuler method), the Jacobian is supposed to be evaluated at ts->ptime like this: 303 if(th->endpoint) { 304 ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 305 } 306 but ts->ptime and ts->vec_sol have the same values as th->stage_time and th->X in this case. So the code is simplified here. 307 */ 308 for (nadj=0; nadj<ts->numberadjs; nadj++) { 309 ierr = MatMultTranspose(J,VecDeltaLam[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 310 ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecSensiTemp[nadj]);CHKERRQ(ierr); 311 } 312 if (ts->vecs_sensip) { 313 ierr = TSRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 314 for (nadj=0; nadj<ts->numberadjs; nadj++) { 315 ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr); 316 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecDeltaMu[nadj]);CHKERRQ(ierr); 317 } 318 } 319 } 320 321 ts->ptime += ts->time_step; 322 ts->steps++; 323 th->status = TS_STEP_COMPLETE; 324 PetscFunctionReturn(0); 325 } 326 327 #undef __FUNCT__ 328 #define __FUNCT__ "TSInterpolate_Theta" 329 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 330 { 331 TS_Theta *th = (TS_Theta*)ts->data; 332 PetscReal alpha = t - ts->ptime; 333 PetscErrorCode ierr; 334 335 PetscFunctionBegin; 336 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 337 if (th->endpoint) alpha *= th->Theta; 338 ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 339 PetscFunctionReturn(0); 340 } 341 342 /*------------------------------------------------------------*/ 343 #undef __FUNCT__ 344 #define __FUNCT__ "TSReset_Theta" 345 static PetscErrorCode TSReset_Theta(TS ts) 346 { 347 TS_Theta *th = (TS_Theta*)ts->data; 348 PetscErrorCode ierr; 349 350 PetscFunctionBegin; 351 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 352 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 353 ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 354 ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 355 if(ts->reverse_mode) { 356 ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr); 357 if(th->VecDeltaMu) { 358 ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr); 359 } 360 ierr = VecDestroyVecs(ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr); 361 } 362 PetscFunctionReturn(0); 363 } 364 365 #undef __FUNCT__ 366 #define __FUNCT__ "TSDestroy_Theta" 367 static PetscErrorCode TSDestroy_Theta(TS ts) 368 { 369 PetscErrorCode ierr; 370 371 PetscFunctionBegin; 372 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 373 ierr = PetscFree(ts->data);CHKERRQ(ierr); 374 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 375 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 376 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 377 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 378 PetscFunctionReturn(0); 379 } 380 381 /* 382 This defines the nonlinear equation that is to be solved with SNES 383 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 384 */ 385 #undef __FUNCT__ 386 #define __FUNCT__ "SNESTSFormFunction_Theta" 387 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 388 { 389 TS_Theta *th = (TS_Theta*)ts->data; 390 PetscErrorCode ierr; 391 Vec X0,Xdot; 392 DM dm,dmsave; 393 PetscReal shift = 1./(th->Theta*ts->time_step); 394 395 PetscFunctionBegin; 396 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 397 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 398 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 399 ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 400 401 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 402 dmsave = ts->dm; 403 ts->dm = dm; 404 ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 405 ts->dm = dmsave; 406 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 407 PetscFunctionReturn(0); 408 } 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "SNESTSFormJacobian_Theta" 412 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 413 { 414 TS_Theta *th = (TS_Theta*)ts->data; 415 PetscErrorCode ierr; 416 Vec Xdot; 417 DM dm,dmsave; 418 PetscReal shift = 1./(th->Theta*ts->time_step); 419 420 PetscFunctionBegin; 421 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 422 423 /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 424 ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 425 426 dmsave = ts->dm; 427 ts->dm = dm; 428 ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 429 ts->dm = dmsave; 430 ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "TSSetUp_Theta" 436 static PetscErrorCode TSSetUp_Theta(TS ts) 437 { 438 TS_Theta *th = (TS_Theta*)ts->data; 439 PetscErrorCode ierr; 440 SNES snes; 441 TSAdapt adapt; 442 DM dm; 443 444 PetscFunctionBegin; 445 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 446 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 447 ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 448 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 449 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 450 if (dm) { 451 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 452 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 453 } 454 if (th->Theta == 0.5 && th->endpoint) th->order = 2; 455 else th->order = 1; 456 457 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 458 if (!th->adapt) { 459 ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr); 460 } 461 if (ts->reverse_mode) { 462 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr); 463 if(ts->vecs_sensip) { 464 ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr); 465 } 466 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr); 467 } 468 PetscFunctionReturn(0); 469 } 470 /*------------------------------------------------------------*/ 471 472 #undef __FUNCT__ 473 #define __FUNCT__ "TSSetFromOptions_Theta" 474 static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts) 475 { 476 TS_Theta *th = (TS_Theta*)ts->data; 477 PetscErrorCode ierr; 478 479 PetscFunctionBegin; 480 ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 481 { 482 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 483 ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 484 ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr); 485 ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr); 486 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 487 } 488 ierr = PetscOptionsTail();CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 492 #undef __FUNCT__ 493 #define __FUNCT__ "TSView_Theta" 494 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 495 { 496 TS_Theta *th = (TS_Theta*)ts->data; 497 PetscBool iascii; 498 PetscErrorCode ierr; 499 500 PetscFunctionBegin; 501 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 502 if (iascii) { 503 ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 504 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 505 } 506 if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 507 PetscFunctionReturn(0); 508 } 509 510 #undef __FUNCT__ 511 #define __FUNCT__ "TSThetaGetTheta_Theta" 512 PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 513 { 514 TS_Theta *th = (TS_Theta*)ts->data; 515 516 PetscFunctionBegin; 517 *theta = th->Theta; 518 PetscFunctionReturn(0); 519 } 520 521 #undef __FUNCT__ 522 #define __FUNCT__ "TSThetaSetTheta_Theta" 523 PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 524 { 525 TS_Theta *th = (TS_Theta*)ts->data; 526 527 PetscFunctionBegin; 528 if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 529 th->Theta = theta; 530 PetscFunctionReturn(0); 531 } 532 533 #undef __FUNCT__ 534 #define __FUNCT__ "TSThetaGetEndpoint_Theta" 535 PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 536 { 537 TS_Theta *th = (TS_Theta*)ts->data; 538 539 PetscFunctionBegin; 540 *endpoint = th->endpoint; 541 PetscFunctionReturn(0); 542 } 543 544 #undef __FUNCT__ 545 #define __FUNCT__ "TSThetaSetEndpoint_Theta" 546 PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 547 { 548 TS_Theta *th = (TS_Theta*)ts->data; 549 550 PetscFunctionBegin; 551 th->endpoint = flg; 552 PetscFunctionReturn(0); 553 } 554 555 #if defined(PETSC_HAVE_COMPLEX) 556 #undef __FUNCT__ 557 #define __FUNCT__ "TSComputeLinearStability_Theta" 558 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 559 { 560 PetscComplex z = xr + xi*PETSC_i,f; 561 TS_Theta *th = (TS_Theta*)ts->data; 562 const PetscReal one = 1.0; 563 564 PetscFunctionBegin; 565 f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 566 *yr = PetscRealPartComplex(f); 567 *yi = PetscImaginaryPartComplex(f); 568 PetscFunctionReturn(0); 569 } 570 #endif 571 572 #undef __FUNCT__ 573 #define __FUNCT__ "TSGetStages_Theta" 574 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 575 { 576 TS_Theta *th = (TS_Theta*)ts->data; 577 578 PetscFunctionBegin; 579 *ns = 1; 580 if(Y) { 581 *Y = &(th->X); 582 } 583 PetscFunctionReturn(0); 584 } 585 586 /* ------------------------------------------------------------ */ 587 /*MC 588 TSTHETA - DAE solver using the implicit Theta method 589 590 Level: beginner 591 592 Options Database: 593 -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 594 -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 595 -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 596 597 Notes: 598 $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 599 $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 600 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 601 602 603 604 This method can be applied to DAE. 605 606 This method is cast as a 1-stage implicit Runge-Kutta method. 607 608 .vb 609 Theta | Theta 610 ------------- 611 | 1 612 .ve 613 614 For the default Theta=0.5, this is also known as the implicit midpoint rule. 615 616 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 617 618 .vb 619 0 | 0 0 620 1 | 1-Theta Theta 621 ------------------- 622 | 1-Theta Theta 623 .ve 624 625 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 626 627 To apply a diagonally implicit RK method to DAE, the stage formula 628 629 $ Y_i = X + h sum_j a_ij Y'_j 630 631 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 632 633 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 634 635 M*/ 636 #undef __FUNCT__ 637 #define __FUNCT__ "TSCreate_Theta" 638 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 639 { 640 TS_Theta *th; 641 PetscErrorCode ierr; 642 643 PetscFunctionBegin; 644 ts->ops->reset = TSReset_Theta; 645 ts->ops->destroy = TSDestroy_Theta; 646 ts->ops->view = TSView_Theta; 647 ts->ops->setup = TSSetUp_Theta; 648 ts->ops->step = TSStep_Theta; 649 ts->ops->interpolate = TSInterpolate_Theta; 650 ts->ops->evaluatestep = TSEvaluateStep_Theta; 651 ts->ops->rollback = TSRollBack_Theta; 652 ts->ops->setfromoptions = TSSetFromOptions_Theta; 653 ts->ops->snesfunction = SNESTSFormFunction_Theta; 654 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 655 #if defined(PETSC_HAVE_COMPLEX) 656 ts->ops->linearstability = TSComputeLinearStability_Theta; 657 #endif 658 ts->ops->getstages = TSGetStages_Theta; 659 ts->ops->stepadj = TSStepAdj_Theta; 660 661 ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 662 ts->data = (void*)th; 663 664 th->extrapolate = PETSC_FALSE; 665 th->Theta = 0.5; 666 th->ccfl = 1.0; 667 th->adapt = PETSC_FALSE; 668 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 669 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 670 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 671 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 672 PetscFunctionReturn(0); 673 } 674 675 #undef __FUNCT__ 676 #define __FUNCT__ "TSThetaGetTheta" 677 /*@ 678 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 679 680 Not Collective 681 682 Input Parameter: 683 . ts - timestepping context 684 685 Output Parameter: 686 . theta - stage abscissa 687 688 Note: 689 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 690 691 Level: Advanced 692 693 .seealso: TSThetaSetTheta() 694 @*/ 695 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 696 { 697 PetscErrorCode ierr; 698 699 PetscFunctionBegin; 700 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 701 PetscValidPointer(theta,2); 702 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 703 PetscFunctionReturn(0); 704 } 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "TSThetaSetTheta" 708 /*@ 709 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 710 711 Not Collective 712 713 Input Parameter: 714 + ts - timestepping context 715 - theta - stage abscissa 716 717 Options Database: 718 . -ts_theta_theta <theta> 719 720 Level: Intermediate 721 722 .seealso: TSThetaGetTheta() 723 @*/ 724 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 725 { 726 PetscErrorCode ierr; 727 728 PetscFunctionBegin; 729 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 730 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 731 PetscFunctionReturn(0); 732 } 733 734 #undef __FUNCT__ 735 #define __FUNCT__ "TSThetaGetEndpoint" 736 /*@ 737 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 738 739 Not Collective 740 741 Input Parameter: 742 . ts - timestepping context 743 744 Output Parameter: 745 . endpoint - PETSC_TRUE when using the endpoint variant 746 747 Level: Advanced 748 749 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 750 @*/ 751 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 752 { 753 PetscErrorCode ierr; 754 755 PetscFunctionBegin; 756 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 757 PetscValidPointer(endpoint,2); 758 ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 759 PetscFunctionReturn(0); 760 } 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "TSThetaSetEndpoint" 764 /*@ 765 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 766 767 Not Collective 768 769 Input Parameter: 770 + ts - timestepping context 771 - flg - PETSC_TRUE to use the endpoint variant 772 773 Options Database: 774 . -ts_theta_endpoint <flg> 775 776 Level: Intermediate 777 778 .seealso: TSTHETA, TSCN 779 @*/ 780 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 781 { 782 PetscErrorCode ierr; 783 784 PetscFunctionBegin; 785 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 786 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 787 PetscFunctionReturn(0); 788 } 789 790 /* 791 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 792 * The creation functions for these specializations are below. 793 */ 794 795 #undef __FUNCT__ 796 #define __FUNCT__ "TSView_BEuler" 797 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 798 { 799 PetscErrorCode ierr; 800 801 PetscFunctionBegin; 802 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 803 PetscFunctionReturn(0); 804 } 805 806 /*MC 807 TSBEULER - ODE solver using the implicit backward Euler method 808 809 Level: beginner 810 811 Notes: 812 TSBEULER is equivalent to TSTHETA with Theta=1.0 813 814 $ -ts_type theta -ts_theta_theta 1. 815 816 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 817 818 M*/ 819 #undef __FUNCT__ 820 #define __FUNCT__ "TSCreate_BEuler" 821 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 822 { 823 PetscErrorCode ierr; 824 825 PetscFunctionBegin; 826 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 827 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 828 ts->ops->view = TSView_BEuler; 829 PetscFunctionReturn(0); 830 } 831 832 #undef __FUNCT__ 833 #define __FUNCT__ "TSView_CN" 834 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 835 { 836 PetscErrorCode ierr; 837 838 PetscFunctionBegin; 839 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842 843 /*MC 844 TSCN - ODE solver using the implicit Crank-Nicolson method. 845 846 Level: beginner 847 848 Notes: 849 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 850 851 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 852 853 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 854 855 M*/ 856 #undef __FUNCT__ 857 #define __FUNCT__ "TSCreate_CN" 858 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 859 { 860 PetscErrorCode ierr; 861 862 PetscFunctionBegin; 863 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 864 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 865 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 866 ts->ops->view = TSView_CN; 867 PetscFunctionReturn(0); 868 } 869