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