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 /* context for time stepping */ 11 PetscReal stage_time; 12 Vec X0,X,Xdot; /* Storage for stages and time derivative */ 13 Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 14 PetscReal Theta; 15 PetscReal ptime; 16 PetscReal time_step; 17 PetscInt order; 18 PetscBool endpoint; 19 PetscBool extrapolate; 20 TSStepStatus status; 21 Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 22 23 /* context for sensitivity analysis */ 24 PetscInt num_tlm; /* Total number of tangent linear equations */ 25 Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 26 Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 27 Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 28 Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 29 Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 30 Mat MatFwdSensip0; /* backup for roll-backs due to events */ 31 Vec VecIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 32 Vec *VecsIntegralSensip0; /* backup for roll-backs due to events */ 33 Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 34 Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 35 Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 36 Vec *VecsAffine; /* Working vectors to store residuals */ 37 /* context for error estimation */ 38 Vec vec_sol_prev; 39 Vec vec_lte_work; 40 } TS_Theta; 41 42 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 43 { 44 TS_Theta *th = (TS_Theta*)ts->data; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 if (X0) { 49 if (dm && dm != ts->dm) { 50 ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 51 } else *X0 = ts->vec_sol; 52 } 53 if (Xdot) { 54 if (dm && dm != ts->dm) { 55 ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 56 } else *Xdot = th->Xdot; 57 } 58 PetscFunctionReturn(0); 59 } 60 61 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 62 { 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 if (X0) { 67 if (dm && dm != ts->dm) { 68 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 69 } 70 } 71 if (Xdot) { 72 if (dm && dm != ts->dm) { 73 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 74 } 75 } 76 PetscFunctionReturn(0); 77 } 78 79 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 80 { 81 PetscFunctionBegin; 82 PetscFunctionReturn(0); 83 } 84 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 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 104 { 105 PetscFunctionBegin; 106 PetscFunctionReturn(0); 107 } 108 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 static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 131 { 132 TS_Theta *th = (TS_Theta*)ts->data; 133 PetscErrorCode ierr; 134 135 PetscFunctionBegin; 136 if (th->endpoint) { 137 /* Evolve ts->vec_costintegral to compute integrals */ 138 if (th->Theta!=1.0) { 139 ierr = TSComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 140 ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 141 } 142 ierr = TSComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 143 ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 144 } else { 145 ierr = TSComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 146 ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 147 } 148 PetscFunctionReturn(0); 149 } 150 151 static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 152 { 153 TS_Theta *th = (TS_Theta*)ts->data; 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 /* backup cost integral */ 158 ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 159 ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 160 PetscFunctionReturn(0); 161 } 162 163 static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 164 { 165 PetscErrorCode ierr; 166 167 PetscFunctionBegin; 168 ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 169 PetscFunctionReturn(0); 170 } 171 172 static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x) 173 { 174 PetscInt nits,lits; 175 PetscErrorCode ierr; 176 177 PetscFunctionBegin; 178 ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 179 ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 180 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 181 ts->snes_its += nits; ts->ksp_its += lits; 182 PetscFunctionReturn(0); 183 } 184 185 static PetscErrorCode TSStep_Theta(TS ts) 186 { 187 TS_Theta *th = (TS_Theta*)ts->data; 188 PetscInt rejections = 0; 189 PetscBool stageok,accept = PETSC_TRUE; 190 PetscReal next_time_step = ts->time_step; 191 PetscErrorCode ierr; 192 193 PetscFunctionBegin; 194 if (!ts->steprollback) { 195 if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 196 ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 197 } 198 199 th->status = TS_STEP_INCOMPLETE; 200 while (!ts->reason && th->status != TS_STEP_COMPLETE) { 201 202 PetscReal shift = 1/(th->Theta*ts->time_step); 203 th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 204 205 ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 206 if (th->extrapolate && !ts->steprestart) { 207 ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 208 } 209 if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 210 if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 211 ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 212 ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 213 ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 214 } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 215 ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 216 } 217 ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 218 ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 219 ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 220 ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 221 if (!stageok) goto reject_step; 222 223 th->status = TS_STEP_PENDING; 224 if (th->endpoint) { 225 ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 226 } else { 227 ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 228 ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 229 } 230 ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 231 th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 232 if (!accept) { 233 ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 234 ts->time_step = next_time_step; 235 goto reject_step; 236 } 237 238 if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 239 th->ptime = ts->ptime; 240 th->time_step = ts->time_step; 241 } 242 243 ts->ptime += ts->time_step; 244 ts->time_step = next_time_step; 245 break; 246 247 reject_step: 248 ts->reject++; accept = PETSC_FALSE; 249 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 250 ts->reason = TS_DIVERGED_STEP_REJECTED; 251 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 252 } 253 } 254 PetscFunctionReturn(0); 255 } 256 257 static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 258 { 259 TS_Theta *th = (TS_Theta*)ts->data; 260 Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 261 Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 262 PetscInt nadj; 263 Mat J,Jp; 264 KSP ksp; 265 PetscReal shift; 266 PetscScalar *xarr; 267 PetscErrorCode ierr; 268 269 PetscFunctionBegin; 270 th->status = TS_STEP_INCOMPLETE; 271 ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 272 ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 273 274 /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 275 th->stage_time = ts->ptime; /* time_step is negative*/ 276 th->ptime = ts->ptime + ts->time_step; 277 th->time_step = -ts->time_step; 278 279 /* Build RHS for first-order adjoint */ 280 if (ts->vec_costintegral) { /* Cost function has an integral term */ 281 ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 282 } 283 for (nadj=0; nadj<ts->numcost; nadj++) { 284 ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 285 ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */ 286 if (ts->vec_costintegral) { 287 ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 288 } 289 } 290 291 /* Build LHS for first-order adjoint */ 292 shift = 1./th->time_step; 293 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 294 ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 295 296 /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 297 for (nadj=0; nadj<ts->numcost; nadj++) { 298 KSPConvergedReason kspreason; 299 ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 300 ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 301 if (kspreason < 0) { 302 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 303 ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 304 } 305 } 306 307 if (ts->vecs_sensi2) { /* U_{n+1} */ 308 /* Get w1 at t_{n+1} from TLM matrix */ 309 ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 310 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 311 /* lambda_s^T F_UU w_1 */ 312 ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 313 if (ts->vecs_fup) { 314 /* lambda_s^T F_UP w_2 */ 315 ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 316 } 317 for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 318 ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 319 ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 320 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 321 if (ts->vecs_fup) { 322 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 323 } 324 } 325 /* Solve stage equation LHS X = RHS for second-order adjoint */ 326 for (nadj=0; nadj<ts->numcost; nadj++) { 327 KSPConvergedReason kspreason; 328 ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 329 ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 330 if (kspreason < 0) { 331 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 332 ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 333 } 334 } 335 } 336 337 /* Update sensitivities, and evaluate integrals if there is any */ 338 shift = 0.0; 339 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 340 ierr = MatScale(J,-1.);CHKERRQ(ierr); 341 for (nadj=0; nadj<ts->numcost; nadj++) { 342 ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 343 ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr); 344 ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 345 if (ts->vecs_sensi2) { 346 ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 347 ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr); 348 ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 349 } 350 } 351 if (ts->vecs_sensip) { 352 ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */ 353 if (ts->vecs_sensi2p) { 354 if (ts->vecs_fpu) { 355 /* lambda_s^T F_PU w_1 */ 356 ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 357 } 358 if (ts->vecs_fpp) { 359 /* lambda_s^T F_PU w_2 */ 360 ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 361 } 362 } 363 for (nadj=0; nadj<ts->numcost; nadj++) { 364 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 365 ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 366 if (ts->vecs_sensi2p) { 367 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 368 ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 369 if (ts->vecs_fpu) { 370 ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 371 } 372 if (ts->vecs_fpp) { 373 ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 374 } 375 } 376 } 377 if (ts->vec_costintegral) { 378 ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 379 for (nadj=0; nadj<ts->numcost; nadj++) { 380 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 381 } 382 } 383 } 384 385 if (ts->vecs_sensi2) { 386 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 387 ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 388 } 389 th->status = TS_STEP_COMPLETE; 390 PetscFunctionReturn(0); 391 } 392 393 static PetscErrorCode TSAdjointStep_Theta(TS ts) 394 { 395 TS_Theta *th = (TS_Theta*)ts->data; 396 Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 397 Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 398 PetscInt nadj; 399 Mat J,Jp; 400 KSP ksp; 401 PetscReal shift; 402 PetscScalar *xarr; 403 PetscErrorCode ierr; 404 405 PetscFunctionBegin; 406 if (th->Theta == 1.) { 407 ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr); 408 PetscFunctionReturn(0); 409 } 410 th->status = TS_STEP_INCOMPLETE; 411 ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 412 ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 413 414 /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 415 th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 416 th->ptime = ts->ptime + ts->time_step; 417 th->time_step = -ts->time_step; 418 419 /* Build RHS for first-order adjoint */ 420 if (ts->vec_costintegral) { /* Cost function has an integral term */ 421 if (th->endpoint) { 422 ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 423 } else { 424 ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 425 } 426 } 427 for (nadj=0; nadj<ts->numcost; nadj++) { 428 ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 429 ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 430 if (ts->vec_costintegral) { 431 ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 432 } 433 } 434 435 /* Build LHS for first-order adjoint */ 436 shift = 1./(th->Theta*th->time_step); 437 if (th->endpoint) { 438 ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 439 } else { 440 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 441 } 442 ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 443 444 /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 445 for (nadj=0; nadj<ts->numcost; nadj++) { 446 KSPConvergedReason kspreason; 447 ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 448 ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 449 if (kspreason < 0) { 450 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 451 ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 452 } 453 } 454 455 /* Second-order adjoint */ 456 if (ts->vecs_sensi2) { /* U_{n+1} */ 457 if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 458 /* Get w1 at t_{n+1} from TLM matrix */ 459 ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 460 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 461 /* lambda_s^T F_UU w_1 */ 462 ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 463 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 464 ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 465 if (ts->vecs_fup) { 466 /* lambda_s^T F_UP w_2 */ 467 ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 468 } 469 for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 470 ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 471 ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 472 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 473 if (ts->vecs_fup) { 474 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 475 } 476 } 477 /* Solve stage equation LHS X = RHS for second-order adjoint */ 478 for (nadj=0; nadj<ts->numcost; nadj++) { 479 KSPConvergedReason kspreason; 480 ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 481 ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 482 if (kspreason < 0) { 483 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 484 ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 485 } 486 } 487 } 488 489 /* Update sensitivities, and evaluate integrals if there is any */ 490 if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 491 shift = 1./((th->Theta-1.)*th->time_step); 492 ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 493 if (ts->vec_costintegral) { /* R_U at t_n */ 494 ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 495 } 496 for (nadj=0; nadj<ts->numcost; nadj++) { 497 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 498 if (ts->vec_costintegral) { 499 ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 500 } 501 ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 502 } 503 504 /* Second-order adjoint */ 505 if (ts->vecs_sensi2) { /* U_n */ 506 /* Get w1 at t_n from TLM matrix */ 507 ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 508 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 509 /* lambda_s^T F_UU w_1 */ 510 ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 511 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 512 ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 513 if (ts->vecs_fup) { 514 /* lambda_s^T F_UU w_2 */ 515 ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 516 } 517 for (nadj=0; nadj<ts->numcost; nadj++) { 518 /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UU w_2 */ 519 ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 520 ierr = VecAXPY(ts->vecs_sensi2[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 521 if (ts->vecs_fup) { 522 ierr = VecAXPY(ts->vecs_sensi2[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 523 } 524 ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 525 } 526 } 527 528 if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 529 /* U_{n+1} */ 530 ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 531 if (ts->vec_costintegral) { 532 ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 533 } 534 for (nadj=0; nadj<ts->numcost; nadj++) { 535 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 536 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 537 } 538 if (ts->vecs_sensip2) { /* second-order */ 539 /* lambda_s^T F_PU w_1 */ 540 ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 541 /* lambda_s^T F_PP w_2 */ 542 ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 543 for (nadj=0; nadj<ts->numcost; nadj++) { 544 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 545 ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 546 if (ts->vecs_fpu) { 547 ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 548 } 549 if (ts->vecs_fpp) { 550 ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 551 } 552 } 553 } 554 555 /* U_s */ 556 ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 557 if (ts->vec_costintegral) { 558 ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 559 } 560 for (nadj=0; nadj<ts->numcost; nadj++) { 561 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 562 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 563 if (ts->vecs_sensip2) { /* second-order */ 564 /* lambda_s^T F_PU w_1 */ 565 ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 566 /* lambda_s^T F_PP w_2 */ 567 ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 568 for (nadj=0; nadj<ts->numcost; nadj++) { 569 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 570 ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 571 if (ts->vecs_fpu) { 572 ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 573 } 574 if (ts->vecs_fpp) { 575 ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 576 } 577 } 578 } 579 } 580 } 581 } else { /* one-stage case */ 582 shift = 0.0; 583 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 584 if (ts->vec_costintegral) { 585 ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 586 } 587 for (nadj=0; nadj<ts->numcost; nadj++) { 588 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 589 ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 590 if (ts->vec_costintegral) { 591 ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 592 } 593 } 594 if (ts->vecs_sensip) { 595 ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 596 for (nadj=0; nadj<ts->numcost; nadj++) { 597 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 598 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 599 } 600 if (ts->vec_costintegral) { 601 ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 602 for (nadj=0; nadj<ts->numcost; nadj++) { 603 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 604 } 605 } 606 } 607 } 608 609 th->status = TS_STEP_COMPLETE; 610 PetscFunctionReturn(0); 611 } 612 613 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 614 { 615 TS_Theta *th = (TS_Theta*)ts->data; 616 PetscReal dt = t - ts->ptime; 617 PetscErrorCode ierr; 618 619 PetscFunctionBegin; 620 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 621 if (th->endpoint) dt *= th->Theta; 622 ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 623 PetscFunctionReturn(0); 624 } 625 626 static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 627 { 628 TS_Theta *th = (TS_Theta*)ts->data; 629 Vec X = ts->vec_sol; /* X = solution */ 630 Vec Y = th->vec_lte_work; /* Y = X + LTE */ 631 PetscReal wltea,wlter; 632 PetscErrorCode ierr; 633 634 PetscFunctionBegin; 635 if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 636 /* Cannot compute LTE in first step or in restart after event */ 637 if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 638 /* Compute LTE using backward differences with non-constant time step */ 639 { 640 PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 641 PetscReal a = 1 + h_prev/h; 642 PetscScalar scal[3]; Vec vecs[3]; 643 scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 644 vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 645 ierr = VecCopy(X,Y);CHKERRQ(ierr); 646 ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 647 ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 648 } 649 if (order) *order = 2; 650 PetscFunctionReturn(0); 651 } 652 653 static PetscErrorCode TSRollBack_Theta(TS ts) 654 { 655 TS_Theta *th = (TS_Theta*)ts->data; 656 PetscInt ncost; 657 PetscErrorCode ierr; 658 659 PetscFunctionBegin; 660 ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 661 if (ts->vec_costintegral && ts->costintegralfwd) { 662 ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 663 } 664 th->status = TS_STEP_INCOMPLETE; 665 if (ts->mat_sensip) { 666 ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 667 } 668 if (ts->vecs_integral_sensip) { 669 for (ncost=0;ncost<ts->numcost;ncost++) { 670 ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 671 } 672 } 673 PetscFunctionReturn(0); 674 } 675 676 static PetscErrorCode TSForwardStep_Theta(TS ts) 677 { 678 TS_Theta *th = (TS_Theta*)ts->data; 679 Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 680 Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 681 PetscInt ncost,ntlm; 682 KSP ksp; 683 Mat J,Jp; 684 PetscReal shift; 685 PetscScalar *barr,*xarr; 686 PetscErrorCode ierr; 687 688 PetscFunctionBegin; 689 ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 690 691 for (ncost=0; ncost<ts->numcost; ncost++) { 692 if (ts->vecs_integral_sensip) { 693 ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 694 } 695 } 696 697 ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 698 ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 699 700 /* Build RHS */ 701 if (th->endpoint) { /* 2-stage method*/ 702 shift = 1./((th->Theta-1.)*th->time_step); 703 ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 704 ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 705 ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 706 707 /* Add the f_p forcing terms */ 708 if (ts->Jacp) { 709 ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 710 ierr = MatAXPY(MatDeltaFwdSensip,(1.-th->Theta)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 711 ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 712 ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 713 } 714 } else { /* 1-stage method */ 715 shift = 0.0; 716 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 717 ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 718 ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 719 720 /* Add the f_p forcing terms */ 721 if (ts->Jacp) { 722 ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 723 ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 724 } 725 } 726 727 /* Build LHS */ 728 shift = 1/(th->Theta*th->time_step); 729 if (th->endpoint) { 730 ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 731 } else { 732 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 733 } 734 ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 735 736 /* 737 Evaluate the first stage of integral gradients with the 2-stage method: 738 drdu|t_n*S(t_n) + drdp|t_n 739 This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 740 */ 741 if (th->endpoint) { /* 2-stage method only */ 742 if (ts->vecs_integral_sensip) { 743 ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 744 if (ts->vecs_drdp) { 745 ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 746 } 747 for (ncost=0; ncost<ts->numcost; ncost++) { 748 ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 749 if (ts->vecs_drdp) { 750 ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 751 } 752 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 753 } 754 } 755 } 756 757 /* Solve the tangent linear equation for forward sensitivities to parameters */ 758 for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 759 KSPConvergedReason kspreason; 760 ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 761 ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 762 if (th->endpoint) { 763 ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 764 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 765 ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 766 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 767 ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 768 } else { 769 ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 770 } 771 ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 772 if (kspreason < 0) { 773 ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 774 ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 775 } 776 ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 777 ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 778 } 779 780 781 /* 782 Evaluate the second stage of integral gradients with the 2-stage method: 783 drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 784 */ 785 if (ts->vecs_integral_sensip) { 786 if (!th->endpoint) { 787 ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 788 ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 789 if (ts->vecs_drdp) { 790 ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 791 } 792 for (ncost=0; ncost<ts->numcost; ncost++) { 793 ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 794 if (ts->vecs_drdp) { 795 ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 796 } 797 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 798 } 799 ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 800 } else { 801 ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 802 if (ts->vecs_drdp) { 803 ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 804 } 805 for (ncost=0; ncost<ts->numcost; ncost++) { 806 ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 807 if (ts->vecs_drdp) { 808 ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 809 } 810 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 811 } 812 } 813 } else { 814 if (!th->endpoint) { 815 ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 816 } 817 } 818 PetscFunctionReturn(0); 819 } 820 821 static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 822 { 823 TS_Theta *th = (TS_Theta*)ts->data; 824 825 PetscFunctionBegin; 826 if (ns) *ns = 1; 827 if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 828 PetscFunctionReturn(0); 829 } 830 831 /*------------------------------------------------------------*/ 832 static PetscErrorCode TSReset_Theta(TS ts) 833 { 834 TS_Theta *th = (TS_Theta*)ts->data; 835 PetscErrorCode ierr; 836 837 PetscFunctionBegin; 838 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 839 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 840 ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 841 ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 842 843 ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 844 ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 845 846 ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 847 if (ts->forward_solve) { 848 if (ts->vecs_integral_sensip) { 849 ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 850 ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 851 } 852 ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 853 ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 854 ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 855 } 856 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 857 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 858 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 859 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 860 ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 861 ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 862 863 PetscFunctionReturn(0); 864 } 865 866 static PetscErrorCode TSAdjointReset_Theta(TS ts) 867 { 868 TS_Theta *th = (TS_Theta*)ts->data; 869 PetscErrorCode ierr; 870 871 PetscFunctionBegin; 872 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 873 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 874 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 875 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 876 ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 877 ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 878 PetscFunctionReturn(0); 879 } 880 881 static PetscErrorCode TSDestroy_Theta(TS ts) 882 { 883 PetscErrorCode ierr; 884 885 PetscFunctionBegin; 886 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 887 if (ts->dm) { 888 ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 889 ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 890 } 891 ierr = PetscFree(ts->data);CHKERRQ(ierr); 892 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 893 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 894 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 895 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 896 PetscFunctionReturn(0); 897 } 898 899 /* 900 This defines the nonlinear equation that is to be solved with SNES 901 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 902 */ 903 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 904 { 905 TS_Theta *th = (TS_Theta*)ts->data; 906 PetscErrorCode ierr; 907 Vec X0,Xdot; 908 DM dm,dmsave; 909 PetscReal shift = 1/(th->Theta*ts->time_step); 910 911 PetscFunctionBegin; 912 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 913 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 914 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 915 ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 916 917 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 918 dmsave = ts->dm; 919 ts->dm = dm; 920 ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 921 ts->dm = dmsave; 922 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 923 PetscFunctionReturn(0); 924 } 925 926 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 927 { 928 TS_Theta *th = (TS_Theta*)ts->data; 929 PetscErrorCode ierr; 930 Vec Xdot; 931 DM dm,dmsave; 932 PetscReal shift = 1/(th->Theta*ts->time_step); 933 934 PetscFunctionBegin; 935 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 936 /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 937 ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 938 939 dmsave = ts->dm; 940 ts->dm = dm; 941 ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 942 ts->dm = dmsave; 943 ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 944 PetscFunctionReturn(0); 945 } 946 947 static PetscErrorCode TSForwardSetUp_Theta(TS ts) 948 { 949 TS_Theta *th = (TS_Theta*)ts->data; 950 PetscErrorCode ierr; 951 952 PetscFunctionBegin; 953 /* combine sensitivities to parameters and sensitivities to initial values into one array */ 954 th->num_tlm = ts->num_parameters; 955 ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 956 if (ts->vecs_integral_sensip) { 957 ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 958 } 959 /* backup sensitivity results for roll-backs */ 960 ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 961 962 if (ts->vecs_integral_sensip) { 963 ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 964 } 965 ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 966 ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 967 PetscFunctionReturn(0); 968 } 969 970 static PetscErrorCode TSForwardReset_Theta(TS ts) 971 { 972 TS_Theta *th = (TS_Theta*)ts->data; 973 PetscErrorCode ierr; 974 975 PetscFunctionBegin; 976 if (ts->vecs_integral_sensip) { 977 ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 978 ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 979 } 980 ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 981 ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 982 ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 983 PetscFunctionReturn(0); 984 } 985 986 static PetscErrorCode TSSetUp_Theta(TS ts) 987 { 988 TS_Theta *th = (TS_Theta*)ts->data; 989 PetscBool match; 990 PetscErrorCode ierr; 991 992 PetscFunctionBegin; 993 if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 994 ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 995 } 996 if (!th->X) { 997 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 998 } 999 if (!th->Xdot) { 1000 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 1001 } 1002 if (!th->X0) { 1003 ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 1004 } 1005 if (th->endpoint) { 1006 ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 1007 } 1008 1009 th->order = (th->Theta == 0.5) ? 2 : 1; 1010 1011 ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 1012 ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 1013 ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 1014 1015 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 1016 ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 1017 ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 1018 if (!match) { 1019 ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 1020 ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 1021 } 1022 ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1023 PetscFunctionReturn(0); 1024 } 1025 1026 /*------------------------------------------------------------*/ 1027 1028 static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1029 { 1030 TS_Theta *th = (TS_Theta*)ts->data; 1031 PetscErrorCode ierr; 1032 1033 PetscFunctionBegin; 1034 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1035 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 1036 if (ts->vecs_sensip) { 1037 ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 1038 } 1039 if (ts->vecs_sensi2) { 1040 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1041 ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 1042 } 1043 if (ts->vecs_sensi2p) { 1044 ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 1045 } 1046 PetscFunctionReturn(0); 1047 } 1048 1049 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1050 { 1051 TS_Theta *th = (TS_Theta*)ts->data; 1052 PetscErrorCode ierr; 1053 1054 PetscFunctionBegin; 1055 ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1056 { 1057 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 1058 ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr); 1059 ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 1060 } 1061 ierr = PetscOptionsTail();CHKERRQ(ierr); 1062 PetscFunctionReturn(0); 1063 } 1064 1065 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1066 { 1067 TS_Theta *th = (TS_Theta*)ts->data; 1068 PetscBool iascii; 1069 PetscErrorCode ierr; 1070 1071 PetscFunctionBegin; 1072 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1073 if (iascii) { 1074 ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1075 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1076 } 1077 PetscFunctionReturn(0); 1078 } 1079 1080 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 1081 { 1082 TS_Theta *th = (TS_Theta*)ts->data; 1083 1084 PetscFunctionBegin; 1085 *theta = th->Theta; 1086 PetscFunctionReturn(0); 1087 } 1088 1089 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 1090 { 1091 TS_Theta *th = (TS_Theta*)ts->data; 1092 1093 PetscFunctionBegin; 1094 if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 1095 th->Theta = theta; 1096 th->order = (th->Theta == 0.5) ? 2 : 1; 1097 PetscFunctionReturn(0); 1098 } 1099 1100 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 1101 { 1102 TS_Theta *th = (TS_Theta*)ts->data; 1103 1104 PetscFunctionBegin; 1105 *endpoint = th->endpoint; 1106 PetscFunctionReturn(0); 1107 } 1108 1109 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1110 { 1111 TS_Theta *th = (TS_Theta*)ts->data; 1112 1113 PetscFunctionBegin; 1114 th->endpoint = flg; 1115 PetscFunctionReturn(0); 1116 } 1117 1118 #if defined(PETSC_HAVE_COMPLEX) 1119 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1120 { 1121 PetscComplex z = xr + xi*PETSC_i,f; 1122 TS_Theta *th = (TS_Theta*)ts->data; 1123 const PetscReal one = 1.0; 1124 1125 PetscFunctionBegin; 1126 f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1127 *yr = PetscRealPartComplex(f); 1128 *yi = PetscImaginaryPartComplex(f); 1129 PetscFunctionReturn(0); 1130 } 1131 #endif 1132 1133 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 1134 { 1135 TS_Theta *th = (TS_Theta*)ts->data; 1136 1137 PetscFunctionBegin; 1138 if (ns) *ns = 1; 1139 if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 1140 PetscFunctionReturn(0); 1141 } 1142 1143 /* ------------------------------------------------------------ */ 1144 /*MC 1145 TSTHETA - DAE solver using the implicit Theta method 1146 1147 Level: beginner 1148 1149 Options Database: 1150 + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1151 . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 1152 - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 1153 1154 Notes: 1155 $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 1156 $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 1157 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 1158 1159 This method can be applied to DAE. 1160 1161 This method is cast as a 1-stage implicit Runge-Kutta method. 1162 1163 .vb 1164 Theta | Theta 1165 ------------- 1166 | 1 1167 .ve 1168 1169 For the default Theta=0.5, this is also known as the implicit midpoint rule. 1170 1171 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1172 1173 .vb 1174 0 | 0 0 1175 1 | 1-Theta Theta 1176 ------------------- 1177 | 1-Theta Theta 1178 .ve 1179 1180 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1181 1182 To apply a diagonally implicit RK method to DAE, the stage formula 1183 1184 $ Y_i = X + h sum_j a_ij Y'_j 1185 1186 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1187 1188 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1189 1190 M*/ 1191 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1192 { 1193 TS_Theta *th; 1194 PetscErrorCode ierr; 1195 1196 PetscFunctionBegin; 1197 ts->ops->reset = TSReset_Theta; 1198 ts->ops->adjointreset = TSAdjointReset_Theta; 1199 ts->ops->destroy = TSDestroy_Theta; 1200 ts->ops->view = TSView_Theta; 1201 ts->ops->setup = TSSetUp_Theta; 1202 ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1203 ts->ops->adjointreset = TSAdjointReset_Theta; 1204 ts->ops->step = TSStep_Theta; 1205 ts->ops->interpolate = TSInterpolate_Theta; 1206 ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 1207 ts->ops->rollback = TSRollBack_Theta; 1208 ts->ops->setfromoptions = TSSetFromOptions_Theta; 1209 ts->ops->snesfunction = SNESTSFormFunction_Theta; 1210 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1211 #if defined(PETSC_HAVE_COMPLEX) 1212 ts->ops->linearstability = TSComputeLinearStability_Theta; 1213 #endif 1214 ts->ops->getstages = TSGetStages_Theta; 1215 ts->ops->adjointstep = TSAdjointStep_Theta; 1216 ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1217 ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 1218 ts->default_adapt_type = TSADAPTNONE; 1219 1220 ts->ops->forwardsetup = TSForwardSetUp_Theta; 1221 ts->ops->forwardreset = TSForwardReset_Theta; 1222 ts->ops->forwardstep = TSForwardStep_Theta; 1223 ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1224 1225 ts->usessnes = PETSC_TRUE; 1226 1227 ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1228 ts->data = (void*)th; 1229 1230 th->VecsDeltaLam = NULL; 1231 th->VecsDeltaMu = NULL; 1232 th->VecsSensiTemp = NULL; 1233 th->VecsSensi2Temp = NULL; 1234 1235 th->extrapolate = PETSC_FALSE; 1236 th->Theta = 0.5; 1237 th->order = 2; 1238 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1239 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1240 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1241 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1242 PetscFunctionReturn(0); 1243 } 1244 1245 /*@ 1246 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 1247 1248 Not Collective 1249 1250 Input Parameter: 1251 . ts - timestepping context 1252 1253 Output Parameter: 1254 . theta - stage abscissa 1255 1256 Note: 1257 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 1258 1259 Level: Advanced 1260 1261 .seealso: TSThetaSetTheta() 1262 @*/ 1263 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 1264 { 1265 PetscErrorCode ierr; 1266 1267 PetscFunctionBegin; 1268 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1269 PetscValidPointer(theta,2); 1270 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 1271 PetscFunctionReturn(0); 1272 } 1273 1274 /*@ 1275 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 1276 1277 Not Collective 1278 1279 Input Parameter: 1280 + ts - timestepping context 1281 - theta - stage abscissa 1282 1283 Options Database: 1284 . -ts_theta_theta <theta> 1285 1286 Level: Intermediate 1287 1288 .seealso: TSThetaGetTheta() 1289 @*/ 1290 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 1291 { 1292 PetscErrorCode ierr; 1293 1294 PetscFunctionBegin; 1295 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1296 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 1297 PetscFunctionReturn(0); 1298 } 1299 1300 /*@ 1301 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1302 1303 Not Collective 1304 1305 Input Parameter: 1306 . ts - timestepping context 1307 1308 Output Parameter: 1309 . endpoint - PETSC_TRUE when using the endpoint variant 1310 1311 Level: Advanced 1312 1313 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 1314 @*/ 1315 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 1316 { 1317 PetscErrorCode ierr; 1318 1319 PetscFunctionBegin; 1320 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1321 PetscValidPointer(endpoint,2); 1322 ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 1323 PetscFunctionReturn(0); 1324 } 1325 1326 /*@ 1327 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1328 1329 Not Collective 1330 1331 Input Parameter: 1332 + ts - timestepping context 1333 - flg - PETSC_TRUE to use the endpoint variant 1334 1335 Options Database: 1336 . -ts_theta_endpoint <flg> 1337 1338 Level: Intermediate 1339 1340 .seealso: TSTHETA, TSCN 1341 @*/ 1342 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1343 { 1344 PetscErrorCode ierr; 1345 1346 PetscFunctionBegin; 1347 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1348 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1349 PetscFunctionReturn(0); 1350 } 1351 1352 /* 1353 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1354 * The creation functions for these specializations are below. 1355 */ 1356 1357 static PetscErrorCode TSSetUp_BEuler(TS ts) 1358 { 1359 TS_Theta *th = (TS_Theta*)ts->data; 1360 PetscErrorCode ierr; 1361 1362 PetscFunctionBegin; 1363 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"); 1364 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"); 1365 ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 1366 PetscFunctionReturn(0); 1367 } 1368 1369 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1370 { 1371 PetscFunctionBegin; 1372 PetscFunctionReturn(0); 1373 } 1374 1375 /*MC 1376 TSBEULER - ODE solver using the implicit backward Euler method 1377 1378 Level: beginner 1379 1380 Notes: 1381 TSBEULER is equivalent to TSTHETA with Theta=1.0 1382 1383 $ -ts_type theta -ts_theta_theta 1.0 1384 1385 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1386 1387 M*/ 1388 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1389 { 1390 PetscErrorCode ierr; 1391 1392 PetscFunctionBegin; 1393 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1394 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 1395 ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 1396 ts->ops->setup = TSSetUp_BEuler; 1397 ts->ops->view = TSView_BEuler; 1398 PetscFunctionReturn(0); 1399 } 1400 1401 static PetscErrorCode TSSetUp_CN(TS ts) 1402 { 1403 TS_Theta *th = (TS_Theta*)ts->data; 1404 PetscErrorCode ierr; 1405 1406 PetscFunctionBegin; 1407 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"); 1408 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"); 1409 ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 1410 PetscFunctionReturn(0); 1411 } 1412 1413 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1414 { 1415 PetscFunctionBegin; 1416 PetscFunctionReturn(0); 1417 } 1418 1419 /*MC 1420 TSCN - ODE solver using the implicit Crank-Nicolson method. 1421 1422 Level: beginner 1423 1424 Notes: 1425 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 1426 1427 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1428 1429 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1430 1431 M*/ 1432 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1433 { 1434 PetscErrorCode ierr; 1435 1436 PetscFunctionBegin; 1437 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1438 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1439 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 1440 ts->ops->setup = TSSetUp_CN; 1441 ts->ops->view = TSView_CN; 1442 PetscFunctionReturn(0); 1443 } 1444