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