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 /* context for time stepping */ 11 PetscReal stage_time; 12 Vec X0,X,Xdot; /* Storage for stages and time derivative */ 13 Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 14 PetscReal Theta; 15 PetscReal ptime; 16 PetscReal time_step; 17 PetscInt order; 18 PetscBool endpoint; 19 PetscBool extrapolate; 20 TSStepStatus status; 21 Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 22 23 /* context for sensitivity analysis */ 24 Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 25 Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 26 Vec *VecsSensiTemp; /* Vector to be timed with Jacobian transpose */ 27 Vec *VecsDeltaFwdSensi; /* Increment of the forward sensitivity at stage */ 28 Vec *VecsDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 29 Vec VecIntegralSensiTemp; /* Working vector for forward integral sensitivity */ 30 Vec VecIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 31 32 /* context for error estimation */ 33 Vec vec_sol_prev; 34 Vec vec_lte_work; 35 } TS_Theta; 36 37 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 38 { 39 TS_Theta *th = (TS_Theta*)ts->data; 40 PetscErrorCode ierr; 41 42 PetscFunctionBegin; 43 if (X0) { 44 if (dm && dm != ts->dm) { 45 ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 46 } else *X0 = ts->vec_sol; 47 } 48 if (Xdot) { 49 if (dm && dm != ts->dm) { 50 ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 51 } else *Xdot = th->Xdot; 52 } 53 PetscFunctionReturn(0); 54 } 55 56 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 57 { 58 PetscErrorCode ierr; 59 60 PetscFunctionBegin; 61 if (X0) { 62 if (dm && dm != ts->dm) { 63 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 64 } 65 } 66 if (Xdot) { 67 if (dm && dm != ts->dm) { 68 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 69 } 70 } 71 PetscFunctionReturn(0); 72 } 73 74 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 75 { 76 PetscFunctionBegin; 77 PetscFunctionReturn(0); 78 } 79 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 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 99 { 100 PetscFunctionBegin; 101 PetscFunctionReturn(0); 102 } 103 104 static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 105 { 106 TS ts = (TS)ctx; 107 PetscErrorCode ierr; 108 Vec X0,Xdot,X0_sub,Xdot_sub; 109 110 PetscFunctionBegin; 111 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 112 ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 113 114 ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 115 ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 116 117 ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 118 ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 119 120 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 121 ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 122 PetscFunctionReturn(0); 123 } 124 125 static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 126 { 127 TS_Theta *th = (TS_Theta*)ts->data; 128 PetscErrorCode ierr; 129 130 PetscFunctionBegin; 131 /* backup cost integral */ 132 ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 133 if (th->endpoint) { 134 ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 135 ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 136 } 137 ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 138 if (th->endpoint) { 139 ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 140 } else { 141 ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 142 } 143 PetscFunctionReturn(0); 144 } 145 146 static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 147 { 148 TS_Theta *th = (TS_Theta*)ts->data; 149 PetscErrorCode ierr; 150 151 PetscFunctionBegin; 152 if (th->endpoint) { 153 /* Evolve ts->vec_costintegral to compute integrals */ 154 ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 155 ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 156 if (th->Theta!=1) { 157 ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 158 ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr); 159 } 160 } else { 161 ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 162 ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 163 } 164 PetscFunctionReturn(0); 165 } 166 167 static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 168 { 169 PetscInt nits,lits; 170 PetscErrorCode ierr; 171 172 PetscFunctionBegin; 173 ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 174 ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 175 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 176 ts->snes_its += nits; ts->ksp_its += lits; 177 PetscFunctionReturn(0); 178 } 179 180 static PetscErrorCode TSStep_Theta(TS ts) 181 { 182 TS_Theta *th = (TS_Theta*)ts->data; 183 PetscInt rejections = 0; 184 PetscBool stageok,accept = PETSC_TRUE; 185 PetscReal next_time_step = ts->time_step; 186 PetscErrorCode ierr; 187 188 PetscFunctionBegin; 189 if (!ts->steprollback) { 190 if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 191 ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 192 } 193 194 th->status = TS_STEP_INCOMPLETE; 195 while (!ts->reason && th->status != TS_STEP_COMPLETE) { 196 197 PetscReal shift = 1/(th->Theta*ts->time_step); 198 th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 199 200 ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 201 if (th->extrapolate && !ts->steprestart) { 202 ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 203 } 204 if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 205 if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 206 ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 207 ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 208 ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 209 } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 210 ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 211 } 212 ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 213 ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 214 ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 215 ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 216 if (!stageok) goto reject_step; 217 218 th->status = TS_STEP_PENDING; 219 if (th->endpoint) { 220 ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 221 } else { 222 ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 223 ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 224 } 225 ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 226 th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 227 if (!accept) { 228 ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 229 ts->time_step = next_time_step; 230 goto reject_step; 231 } 232 233 if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 234 th->ptime = ts->ptime; 235 th->time_step = ts->time_step; 236 } 237 238 ts->ptime += ts->time_step; 239 ts->time_step = next_time_step; 240 break; 241 242 reject_step: 243 ts->reject++; accept = PETSC_FALSE; 244 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 245 ts->reason = TS_DIVERGED_STEP_REJECTED; 246 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 247 } 248 } 249 PetscFunctionReturn(0); 250 } 251 252 static PetscErrorCode TSAdjointStep_Theta(TS ts) 253 { 254 TS_Theta *th = (TS_Theta*)ts->data; 255 Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 256 PetscInt nadj; 257 PetscErrorCode ierr; 258 Mat J,Jp; 259 KSP ksp; 260 PetscReal shift; 261 262 PetscFunctionBegin; 263 th->status = TS_STEP_INCOMPLETE; 264 ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 265 ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 266 267 /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 268 th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/ 269 th->ptime = ts->ptime + ts->time_step; 270 271 /* Build RHS */ 272 if (ts->vec_costintegral) { /* Cost function has an integral term */ 273 if (th->endpoint) { 274 ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr); 275 } else { 276 ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 277 } 278 } 279 for (nadj=0; nadj<ts->numcost; nadj++) { 280 ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 281 ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr); 282 if (ts->vec_costintegral) { 283 ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 284 } 285 } 286 287 /* Build LHS */ 288 shift = -1./(th->Theta*ts->time_step); 289 if (th->endpoint) { 290 ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 291 } else { 292 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 293 } 294 ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 295 296 /* Solve LHS X = RHS */ 297 for (nadj=0; nadj<ts->numcost; nadj++) { 298 ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 299 } 300 301 /* Update sensitivities, and evaluate integrals if there is any */ 302 if(th->endpoint) { /* two-stage case */ 303 if (th->Theta!=1.) { 304 shift = -1./((th->Theta-1.)*ts->time_step); 305 ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 306 if (ts->vec_costintegral) { 307 ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr); 308 } 309 for (nadj=0; nadj<ts->numcost; nadj++) { 310 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 311 if (ts->vec_costintegral) { 312 ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 313 } 314 ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 315 } 316 } else { /* backward Euler */ 317 shift = 0.0; 318 ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 319 for (nadj=0; nadj<ts->numcost; nadj++) { 320 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 321 ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 322 if (ts->vec_costintegral) { 323 ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 324 } 325 } 326 } 327 328 if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 329 ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 330 for (nadj=0; nadj<ts->numcost; nadj++) { 331 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 332 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 333 } 334 if (th->Theta!=1.) { 335 ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 336 for (nadj=0; nadj<ts->numcost; nadj++) { 337 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 338 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 339 } 340 } 341 if (ts->vec_costintegral) { 342 ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 343 for (nadj=0; nadj<ts->numcost; nadj++) { 344 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 345 } 346 if (th->Theta!=1.) { 347 ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 348 for (nadj=0; nadj<ts->numcost; nadj++) { 349 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 350 } 351 } 352 } 353 } 354 } else { /* one-stage case */ 355 shift = 0.0; 356 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 357 if (ts->vec_costintegral) { 358 ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 359 } 360 for (nadj=0; nadj<ts->numcost; nadj++) { 361 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 362 ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 363 if (ts->vec_costintegral) { 364 ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 365 } 366 } 367 if (ts->vecs_sensip) { 368 ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 369 for (nadj=0; nadj<ts->numcost; nadj++) { 370 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 371 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 372 } 373 if (ts->vec_costintegral) { 374 ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 375 for (nadj=0; nadj<ts->numcost; nadj++) { 376 ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 377 } 378 } 379 } 380 } 381 382 th->status = TS_STEP_COMPLETE; 383 PetscFunctionReturn(0); 384 } 385 386 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 387 { 388 TS_Theta *th = (TS_Theta*)ts->data; 389 PetscReal dt = t - ts->ptime; 390 PetscErrorCode ierr; 391 392 PetscFunctionBegin; 393 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 394 if (th->endpoint) dt *= th->Theta; 395 ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 396 PetscFunctionReturn(0); 397 } 398 399 static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 400 { 401 TS_Theta *th = (TS_Theta*)ts->data; 402 Vec X = ts->vec_sol; /* X = solution */ 403 Vec Y = th->vec_lte_work; /* Y = X + LTE */ 404 PetscReal wltea,wlter; 405 PetscErrorCode ierr; 406 407 PetscFunctionBegin; 408 if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 409 /* Cannot compute LTE in first step or in restart after event */ 410 if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 411 /* Compute LTE using backward differences with non-constant time step */ 412 { 413 PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 414 PetscReal a = 1 + h_prev/h; 415 PetscScalar scal[3]; Vec vecs[3]; 416 scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 417 vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 418 ierr = VecCopy(X,Y);CHKERRQ(ierr); 419 ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 420 ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 421 } 422 if (order) *order = 2; 423 PetscFunctionReturn(0); 424 } 425 426 static PetscErrorCode TSRollBack_Theta(TS ts) 427 { 428 TS_Theta *th = (TS_Theta*)ts->data; 429 PetscErrorCode ierr; 430 431 PetscFunctionBegin; 432 ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 433 if (ts->vec_costintegral && ts->costintegralfwd) { 434 ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 435 } 436 th->status = TS_STEP_INCOMPLETE; 437 PetscFunctionReturn(0); 438 } 439 440 static PetscErrorCode TSThetaIntegrandDerivative(TS ts,PetscInt nump,Vec DRncostDY,Vec* s,Vec DRncostDP,Vec VecIntegrandDerivative) 441 { 442 PetscInt ndir; 443 PetscScalar *v; 444 PetscErrorCode ierr; 445 446 PetscFunctionBegin; 447 448 ierr = VecGetArray(VecIntegrandDerivative,&v);CHKERRQ(ierr); 449 for (ndir=0; ndir<nump; ndir++) { 450 ierr = VecDot(DRncostDY,s[ndir],&v[ndir]);CHKERRQ(ierr); 451 } 452 ierr = VecRestoreArray(VecIntegrandDerivative,&v);CHKERRQ(ierr); 453 if (DRncostDP) { 454 ierr = VecAXPY(VecIntegrandDerivative,1.,DRncostDP);CHKERRQ(ierr); 455 } 456 PetscFunctionReturn(0); 457 } 458 459 static PetscErrorCode TSForwardStep_Theta(TS ts) 460 { 461 TS_Theta *th = (TS_Theta*)ts->data; 462 Vec *VecsDeltaFwdSensi = th->VecsDeltaFwdSensi; 463 Vec *VecsDeltaFwdSensip = th->VecsDeltaFwdSensip; 464 PetscInt ndir,ncost; 465 KSP ksp; 466 Mat J,Jp; 467 PetscReal shift; 468 PetscErrorCode ierr; 469 470 PetscFunctionBegin; 471 ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 472 ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 473 474 /* Build RHS */ 475 if (th->endpoint) { /* 2-stage method*/ 476 shift = 1./((th->Theta-1.)*th->time_step); 477 ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 478 for (ndir=0; ndir<ts->num_parameters; ndir++) { 479 ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]); 480 ierr = VecScale(VecsDeltaFwdSensip[ndir],(th->Theta-1.)/th->Theta); 481 } 482 for (ndir=0; ndir<ts->num_initialvalues; ndir++) { 483 ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]); 484 ierr = VecScale(VecsDeltaFwdSensi[ndir],(th->Theta-1.)/th->Theta); 485 } 486 487 if (ts->num_parameters) { /* Add the f_p forcing terms */ 488 ierr = TSForwardComputeRHSJacobianP(ts,th->ptime,th->X0,ts->vecs_jacp);CHKERRQ(ierr); 489 for (ndir=0; ndir<ts->num_parameters; ndir++) { 490 ierr = VecAXPY(VecsDeltaFwdSensip[ndir],(1.-th->Theta)/th->Theta,ts->vecs_jacp[ndir]);CHKERRQ(ierr); 491 } 492 ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->vecs_jacp);CHKERRQ(ierr); 493 for (ndir=0; ndir<ts->num_parameters; ndir++) { 494 ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr); 495 } 496 } 497 } else { /* 1-stage method */ 498 shift = 0.0; 499 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 500 for (ndir=0; ndir<ts->num_parameters; ndir++) { 501 ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]); 502 ierr = VecScale(VecsDeltaFwdSensip[ndir],-1); 503 } 504 for (ndir=0; ndir<ts->num_initialvalues; ndir++) { 505 ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]); 506 ierr = VecScale(VecsDeltaFwdSensi[ndir],-1); 507 } 508 if (ts->num_parameters) { /* Add the f_p forcing terms */ 509 ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,th->X,ts->vecs_jacp);CHKERRQ(ierr); 510 for (ndir=0; ndir<ts->num_parameters; ndir++) { 511 ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr); 512 } 513 } 514 } 515 516 /* Build LHS */ 517 shift = 1/(th->Theta*th->time_step); 518 if (th->endpoint) { 519 ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 520 } else { 521 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 522 } 523 ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 524 525 /* 526 Evaluate the first stage of integral gradients with the 2-stage method: 527 drdy|t_n*S(t_n) + drdp|t_n 528 This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 529 It is important to note that sensitivitesi to parameters (sensip) and sensitivities to initial values (sensi) are independent of each other, but integral sensip relies on sensip and integral sensi requires sensi 530 */ 531 if (th->endpoint) { /* 2-stage method only */ 532 if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) { 533 ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr); 534 } 535 if (ts->vecs_integral_sensip) { 536 ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 537 for (ncost=0; ncost<ts->numcost; ncost++) { 538 ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp); 539 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 540 } 541 } 542 if (ts->vecs_integral_sensi) { 543 for (ncost=0; ncost<ts->numcost; ncost++) { 544 ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp); 545 ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensiTemp);CHKERRQ(ierr); 546 } 547 } 548 } 549 550 /* Solve the tangent linear equation for forward sensitivities to parameters */ 551 for (ndir=0; ndir<ts->num_parameters; ndir++) { 552 if (th->endpoint) { 553 ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],ts->vecs_fwdsensip[ndir]); 554 } else { 555 ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],VecsDeltaFwdSensip[ndir]); 556 ierr = VecAXPY(ts->vecs_fwdsensip[ndir],1.,VecsDeltaFwdSensip[ndir]); 557 } 558 } 559 /* Solve the linear system for trajectory sensitivities to initial values */ 560 for (ndir=0; ndir<ts->num_initialvalues; ndir++) { 561 if (th->endpoint) { 562 ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],ts->vecs_fwdsensi[ndir]); 563 } else { 564 ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],VecsDeltaFwdSensi[ndir]); 565 ierr = VecAXPY(ts->vecs_fwdsensi[ndir],1.,VecsDeltaFwdSensi[ndir]); 566 } 567 } 568 569 /*Evaluate the second stage of integral gradients with the 2-stage method: 570 drdy|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 571 */ 572 if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) { 573 ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 574 } 575 if (ts->vecs_integral_sensip) { 576 ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 577 for (ncost=0; ncost<ts->numcost; ncost++) { 578 ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp); 579 if (th->endpoint) { 580 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 581 } else { 582 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 583 } 584 } 585 } 586 if (ts->vecs_integral_sensi) { 587 for (ncost=0; ncost<ts->numcost; ncost++) { 588 ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp); 589 if (th->endpoint) { 590 ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*th->Theta,th->VecIntegralSensiTemp);CHKERRQ(ierr); 591 } else { 592 ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step,th->VecIntegralSensiTemp);CHKERRQ(ierr); 593 } 594 } 595 } 596 597 if (!th->endpoint) { /* VecsDeltaFwdSensip are the increment for the stage values for the 1-stage method */ 598 for (ndir=0; ndir<ts->num_parameters; ndir++) { 599 ierr = VecAXPY(ts->vecs_fwdsensip[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensip[ndir]); 600 } 601 for (ndir=0; ndir<ts->num_initialvalues; ndir++) { 602 ierr = VecAXPY(ts->vecs_fwdsensi[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensi[ndir]); 603 } 604 } 605 606 PetscFunctionReturn(0); 607 } 608 609 /*------------------------------------------------------------*/ 610 static PetscErrorCode TSReset_Theta(TS ts) 611 { 612 TS_Theta *th = (TS_Theta*)ts->data; 613 PetscErrorCode ierr; 614 615 PetscFunctionBegin; 616 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 617 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 618 ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 619 ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 620 621 ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 622 ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 623 624 ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 625 if (ts->forward_solve) { 626 if (ts->vecs_fwdsensi) { 627 ierr = VecDestroyVecs(ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr); } 628 if (ts->vecs_fwdsensip) { 629 ierr = VecDestroyVecs(ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr); 630 } 631 if (ts->vecs_integral_sensi) { 632 ierr = VecDestroy(&th->VecIntegralSensiTemp);CHKERRQ(ierr); 633 } 634 if (ts->vecs_integral_sensip) { 635 ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 636 } 637 } 638 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 639 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 640 ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 641 642 PetscFunctionReturn(0); 643 } 644 645 static PetscErrorCode TSDestroy_Theta(TS ts) 646 { 647 PetscErrorCode ierr; 648 649 PetscFunctionBegin; 650 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 651 ierr = PetscFree(ts->data);CHKERRQ(ierr); 652 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 653 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 654 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 655 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 656 PetscFunctionReturn(0); 657 } 658 659 /* 660 This defines the nonlinear equation that is to be solved with SNES 661 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 662 */ 663 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 664 { 665 TS_Theta *th = (TS_Theta*)ts->data; 666 PetscErrorCode ierr; 667 Vec X0,Xdot; 668 DM dm,dmsave; 669 PetscReal shift = 1/(th->Theta*ts->time_step); 670 671 PetscFunctionBegin; 672 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 673 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 674 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 675 ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 676 677 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 678 dmsave = ts->dm; 679 ts->dm = dm; 680 ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 681 ts->dm = dmsave; 682 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 683 PetscFunctionReturn(0); 684 } 685 686 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 687 { 688 TS_Theta *th = (TS_Theta*)ts->data; 689 PetscErrorCode ierr; 690 Vec Xdot; 691 DM dm,dmsave; 692 PetscReal shift = 1/(th->Theta*ts->time_step); 693 694 PetscFunctionBegin; 695 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 696 /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 697 ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 698 699 dmsave = ts->dm; 700 ts->dm = dm; 701 ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 702 ts->dm = dmsave; 703 ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 704 PetscFunctionReturn(0); 705 } 706 707 static PetscErrorCode TSForwardSetUp_Theta(TS ts) 708 { 709 TS_Theta *th = (TS_Theta*)ts->data; 710 PetscErrorCode ierr; 711 712 PetscFunctionBegin; 713 if (ts->vecs_fwdsensi) { 714 ierr = VecDuplicateVecs(ts->vecs_fwdsensi[0],ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr); 715 } 716 if (ts->vecs_fwdsensip) { 717 ierr = VecDuplicateVecs(ts->vecs_fwdsensip[0],ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr); 718 } 719 if (ts->vecs_integral_sensi) { 720 ierr = VecDuplicate(ts->vecs_integral_sensi[0],&th->VecIntegralSensiTemp);CHKERRQ(ierr); 721 } 722 if (ts->vecs_integral_sensip) { 723 ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 724 } 725 PetscFunctionReturn(0); 726 } 727 728 static PetscErrorCode TSSetUp_Theta(TS ts) 729 { 730 TS_Theta *th = (TS_Theta*)ts->data; 731 PetscBool match; 732 PetscErrorCode ierr; 733 734 PetscFunctionBegin; 735 if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 736 ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 737 } 738 if (!th->X) { 739 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 740 } 741 if (!th->Xdot) { 742 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 743 } 744 if (!th->X0) { 745 ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 746 } 747 if (th->endpoint) { 748 ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 749 } 750 751 th->order = (th->Theta == 0.5) ? 2 : 1; 752 753 ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 754 ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 755 ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 756 757 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 758 ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 759 ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 760 if (!match) { 761 ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 762 ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 763 } 764 ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 765 PetscFunctionReturn(0); 766 } 767 768 /*------------------------------------------------------------*/ 769 770 static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 771 { 772 TS_Theta *th = (TS_Theta*)ts->data; 773 PetscErrorCode ierr; 774 775 PetscFunctionBegin; 776 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 777 if(ts->vecs_sensip) { 778 ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 779 } 780 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 781 PetscFunctionReturn(0); 782 } 783 784 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 785 { 786 TS_Theta *th = (TS_Theta*)ts->data; 787 PetscErrorCode ierr; 788 789 PetscFunctionBegin; 790 ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 791 { 792 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 793 ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr); 794 ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 795 } 796 ierr = PetscOptionsTail();CHKERRQ(ierr); 797 PetscFunctionReturn(0); 798 } 799 800 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 801 { 802 TS_Theta *th = (TS_Theta*)ts->data; 803 PetscBool iascii; 804 PetscErrorCode ierr; 805 806 PetscFunctionBegin; 807 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 808 if (iascii) { 809 ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 810 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 811 } 812 PetscFunctionReturn(0); 813 } 814 815 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 816 { 817 TS_Theta *th = (TS_Theta*)ts->data; 818 819 PetscFunctionBegin; 820 *theta = th->Theta; 821 PetscFunctionReturn(0); 822 } 823 824 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 825 { 826 TS_Theta *th = (TS_Theta*)ts->data; 827 828 PetscFunctionBegin; 829 if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 830 th->Theta = theta; 831 th->order = (th->Theta == 0.5) ? 2 : 1; 832 PetscFunctionReturn(0); 833 } 834 835 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 836 { 837 TS_Theta *th = (TS_Theta*)ts->data; 838 839 PetscFunctionBegin; 840 *endpoint = th->endpoint; 841 PetscFunctionReturn(0); 842 } 843 844 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 845 { 846 TS_Theta *th = (TS_Theta*)ts->data; 847 848 PetscFunctionBegin; 849 th->endpoint = flg; 850 PetscFunctionReturn(0); 851 } 852 853 #if defined(PETSC_HAVE_COMPLEX) 854 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 855 { 856 PetscComplex z = xr + xi*PETSC_i,f; 857 TS_Theta *th = (TS_Theta*)ts->data; 858 const PetscReal one = 1.0; 859 860 PetscFunctionBegin; 861 f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 862 *yr = PetscRealPartComplex(f); 863 *yi = PetscImaginaryPartComplex(f); 864 PetscFunctionReturn(0); 865 } 866 #endif 867 868 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 869 { 870 TS_Theta *th = (TS_Theta*)ts->data; 871 872 PetscFunctionBegin; 873 if (ns) *ns = 1; 874 if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 875 PetscFunctionReturn(0); 876 } 877 878 /* ------------------------------------------------------------ */ 879 /*MC 880 TSTHETA - DAE solver using the implicit Theta method 881 882 Level: beginner 883 884 Options Database: 885 + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 886 . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 887 - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 888 889 Notes: 890 $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 891 $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 892 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 893 894 This method can be applied to DAE. 895 896 This method is cast as a 1-stage implicit Runge-Kutta method. 897 898 .vb 899 Theta | Theta 900 ------------- 901 | 1 902 .ve 903 904 For the default Theta=0.5, this is also known as the implicit midpoint rule. 905 906 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 907 908 .vb 909 0 | 0 0 910 1 | 1-Theta Theta 911 ------------------- 912 | 1-Theta Theta 913 .ve 914 915 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 916 917 To apply a diagonally implicit RK method to DAE, the stage formula 918 919 $ Y_i = X + h sum_j a_ij Y'_j 920 921 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 922 923 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 924 925 M*/ 926 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 927 { 928 TS_Theta *th; 929 PetscErrorCode ierr; 930 931 PetscFunctionBegin; 932 ts->ops->reset = TSReset_Theta; 933 ts->ops->destroy = TSDestroy_Theta; 934 ts->ops->view = TSView_Theta; 935 ts->ops->setup = TSSetUp_Theta; 936 ts->ops->adjointsetup = TSAdjointSetUp_Theta; 937 ts->ops->step = TSStep_Theta; 938 ts->ops->interpolate = TSInterpolate_Theta; 939 ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 940 ts->ops->rollback = TSRollBack_Theta; 941 ts->ops->setfromoptions = TSSetFromOptions_Theta; 942 ts->ops->snesfunction = SNESTSFormFunction_Theta; 943 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 944 #if defined(PETSC_HAVE_COMPLEX) 945 ts->ops->linearstability = TSComputeLinearStability_Theta; 946 #endif 947 ts->ops->getstages = TSGetStages_Theta; 948 ts->ops->adjointstep = TSAdjointStep_Theta; 949 ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 950 ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 951 ts->default_adapt_type = TSADAPTNONE; 952 ts->ops->forwardsetup = TSForwardSetUp_Theta; 953 ts->ops->forwardstep = TSForwardStep_Theta; 954 955 ts->usessnes = PETSC_TRUE; 956 957 ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 958 ts->data = (void*)th; 959 960 th->VecsDeltaLam = NULL; 961 th->VecsDeltaMu = NULL; 962 th->VecsSensiTemp = NULL; 963 964 th->extrapolate = PETSC_FALSE; 965 th->Theta = 0.5; 966 th->order = 2; 967 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 968 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 969 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 970 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 971 PetscFunctionReturn(0); 972 } 973 974 /*@ 975 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 976 977 Not Collective 978 979 Input Parameter: 980 . ts - timestepping context 981 982 Output Parameter: 983 . theta - stage abscissa 984 985 Note: 986 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 987 988 Level: Advanced 989 990 .seealso: TSThetaSetTheta() 991 @*/ 992 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 993 { 994 PetscErrorCode ierr; 995 996 PetscFunctionBegin; 997 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 998 PetscValidPointer(theta,2); 999 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 1000 PetscFunctionReturn(0); 1001 } 1002 1003 /*@ 1004 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 1005 1006 Not Collective 1007 1008 Input Parameter: 1009 + ts - timestepping context 1010 - theta - stage abscissa 1011 1012 Options Database: 1013 . -ts_theta_theta <theta> 1014 1015 Level: Intermediate 1016 1017 .seealso: TSThetaGetTheta() 1018 @*/ 1019 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 1020 { 1021 PetscErrorCode ierr; 1022 1023 PetscFunctionBegin; 1024 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1025 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 1026 PetscFunctionReturn(0); 1027 } 1028 1029 /*@ 1030 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1031 1032 Not Collective 1033 1034 Input Parameter: 1035 . ts - timestepping context 1036 1037 Output Parameter: 1038 . endpoint - PETSC_TRUE when using the endpoint variant 1039 1040 Level: Advanced 1041 1042 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 1043 @*/ 1044 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 1045 { 1046 PetscErrorCode ierr; 1047 1048 PetscFunctionBegin; 1049 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1050 PetscValidPointer(endpoint,2); 1051 ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 1052 PetscFunctionReturn(0); 1053 } 1054 1055 /*@ 1056 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1057 1058 Not Collective 1059 1060 Input Parameter: 1061 + ts - timestepping context 1062 - flg - PETSC_TRUE to use the endpoint variant 1063 1064 Options Database: 1065 . -ts_theta_endpoint <flg> 1066 1067 Level: Intermediate 1068 1069 .seealso: TSTHETA, TSCN 1070 @*/ 1071 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1072 { 1073 PetscErrorCode ierr; 1074 1075 PetscFunctionBegin; 1076 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1077 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1078 PetscFunctionReturn(0); 1079 } 1080 1081 /* 1082 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1083 * The creation functions for these specializations are below. 1084 */ 1085 1086 static PetscErrorCode TSSetUp_BEuler(TS ts) 1087 { 1088 TS_Theta *th = (TS_Theta*)ts->data; 1089 PetscErrorCode ierr; 1090 1091 PetscFunctionBegin; 1092 if (th->Theta != 1.0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler\n"); 1093 if (th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the endpoint form of the Theta methods when using backward Euler\n"); 1094 ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 1095 PetscFunctionReturn(0); 1096 } 1097 1098 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1099 { 1100 PetscFunctionBegin; 1101 PetscFunctionReturn(0); 1102 } 1103 1104 /*MC 1105 TSBEULER - ODE solver using the implicit backward Euler method 1106 1107 Level: beginner 1108 1109 Notes: 1110 TSBEULER is equivalent to TSTHETA with Theta=1.0 1111 1112 $ -ts_type theta -ts_theta_theta 1.0 1113 1114 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1115 1116 M*/ 1117 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1118 { 1119 PetscErrorCode ierr; 1120 1121 PetscFunctionBegin; 1122 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1123 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 1124 ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 1125 ts->ops->setup = TSSetUp_BEuler; 1126 ts->ops->view = TSView_BEuler; 1127 PetscFunctionReturn(0); 1128 } 1129 1130 static PetscErrorCode TSSetUp_CN(TS ts) 1131 { 1132 TS_Theta *th = (TS_Theta*)ts->data; 1133 PetscErrorCode ierr; 1134 1135 PetscFunctionBegin; 1136 if (th->Theta != 0.5) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (0.5) of theta when using Crank-Nicolson\n"); 1137 if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the midpoint form of the Theta methods when using Crank-Nicolson\n"); 1138 ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 1139 PetscFunctionReturn(0); 1140 } 1141 1142 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1143 { 1144 PetscFunctionBegin; 1145 PetscFunctionReturn(0); 1146 } 1147 1148 /*MC 1149 TSCN - ODE solver using the implicit Crank-Nicolson method. 1150 1151 Level: beginner 1152 1153 Notes: 1154 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 1155 1156 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1157 1158 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1159 1160 M*/ 1161 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1162 { 1163 PetscErrorCode ierr; 1164 1165 PetscFunctionBegin; 1166 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1167 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1168 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 1169 ts->ops->setup = TSSetUp_CN; 1170 ts->ops->view = TSView_CN; 1171 PetscFunctionReturn(0); 1172 } 1173