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