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