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