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