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