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