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