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