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