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