1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petscdraw.h> 3 4 PetscLogEvent TS_AdjointStep, TS_ForwardStep, TS_JacobianPEval; 5 6 /* #define TSADJOINT_STAGE */ 7 8 /* ------------------------ Sensitivity Context ---------------------------*/ 9 10 /*@C 11 TSSetRHSJacobianP - Sets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix. 12 13 Logically Collective on TS 14 15 Input Parameters: 16 + ts - TS context obtained from TSCreate() 17 . Amat - JacobianP matrix 18 . func - function 19 - ctx - [optional] user-defined function context 20 21 Calling sequence of func: 22 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 23 + t - current timestep 24 . U - input vector (current ODE solution) 25 . A - output matrix 26 - ctx - [optional] user-defined function context 27 28 Level: intermediate 29 30 Notes: 31 Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 32 33 .seealso: `TSGetRHSJacobianP()` 34 @*/ 35 PetscErrorCode TSSetRHSJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) 36 { 37 PetscFunctionBegin; 38 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 39 PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 40 41 ts->rhsjacobianp = func; 42 ts->rhsjacobianpctx = ctx; 43 if (Amat) { 44 PetscCall(PetscObjectReference((PetscObject)Amat)); 45 PetscCall(MatDestroy(&ts->Jacprhs)); 46 ts->Jacprhs = Amat; 47 } 48 PetscFunctionReturn(0); 49 } 50 51 /*@C 52 TSGetRHSJacobianP - Gets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix. 53 54 Logically Collective on TS 55 56 Input Parameter: 57 . ts - TS context obtained from TSCreate() 58 59 Output Parameters: 60 + Amat - JacobianP matrix 61 . func - function 62 - ctx - [optional] user-defined function context 63 64 Calling sequence of func: 65 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 66 + t - current timestep 67 . U - input vector (current ODE solution) 68 . A - output matrix 69 - ctx - [optional] user-defined function context 70 71 Level: intermediate 72 73 Notes: 74 Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 75 76 .seealso: `TSSetRHSJacobianP()` 77 @*/ 78 PetscErrorCode TSGetRHSJacobianP(TS ts, Mat *Amat, PetscErrorCode (**func)(TS, PetscReal, Vec, Mat, void *), void **ctx) 79 { 80 PetscFunctionBegin; 81 if (func) *func = ts->rhsjacobianp; 82 if (ctx) *ctx = ts->rhsjacobianpctx; 83 if (Amat) *Amat = ts->Jacprhs; 84 PetscFunctionReturn(0); 85 } 86 87 /*@C 88 TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 89 90 Collective on TS 91 92 Input Parameters: 93 . ts - The TS context obtained from TSCreate() 94 95 Level: developer 96 97 .seealso: `TSSetRHSJacobianP()` 98 @*/ 99 PetscErrorCode TSComputeRHSJacobianP(TS ts, PetscReal t, Vec U, Mat Amat) 100 { 101 PetscFunctionBegin; 102 if (!Amat) PetscFunctionReturn(0); 103 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 104 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 105 106 PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 107 PetscFunctionReturn(0); 108 } 109 110 /*@C 111 TSSetIJacobianP - Sets the function that computes the Jacobian of F w.r.t. the parameters P where F(Udot,U,t) = G(U,P,t), as well as the location to store the matrix. 112 113 Logically Collective on TS 114 115 Input Parameters: 116 + ts - TS context obtained from TSCreate() 117 . Amat - JacobianP matrix 118 . func - function 119 - ctx - [optional] user-defined function context 120 121 Calling sequence of func: 122 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 123 + t - current timestep 124 . U - input vector (current ODE solution) 125 . Udot - time derivative of state vector 126 . shift - shift to apply, see note below 127 . A - output matrix 128 - ctx - [optional] user-defined function context 129 130 Level: intermediate 131 132 Notes: 133 Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 134 135 .seealso: 136 @*/ 137 PetscErrorCode TSSetIJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *ctx) 138 { 139 PetscFunctionBegin; 140 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 141 PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 142 143 ts->ijacobianp = func; 144 ts->ijacobianpctx = ctx; 145 if (Amat) { 146 PetscCall(PetscObjectReference((PetscObject)Amat)); 147 PetscCall(MatDestroy(&ts->Jacp)); 148 ts->Jacp = Amat; 149 } 150 PetscFunctionReturn(0); 151 } 152 153 /*@C 154 TSComputeIJacobianP - Runs the user-defined IJacobianP function. 155 156 Collective on TS 157 158 Input Parameters: 159 + ts - the TS context 160 . t - current timestep 161 . U - state vector 162 . Udot - time derivative of state vector 163 . shift - shift to apply, see note below 164 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 165 166 Output Parameters: 167 . A - Jacobian matrix 168 169 Level: developer 170 171 .seealso: `TSSetIJacobianP()` 172 @*/ 173 PetscErrorCode TSComputeIJacobianP(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat Amat, PetscBool imex) 174 { 175 PetscFunctionBegin; 176 if (!Amat) PetscFunctionReturn(0); 177 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 178 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 179 PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4); 180 181 PetscCall(PetscLogEventBegin(TS_JacobianPEval, ts, U, Amat, 0)); 182 if (ts->ijacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->ijacobianp)(ts, t, U, Udot, shift, Amat, ts->ijacobianpctx)); 183 if (imex) { 184 if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */ 185 PetscBool assembled; 186 PetscCall(MatZeroEntries(Amat)); 187 PetscCall(MatAssembled(Amat, &assembled)); 188 if (!assembled) { 189 PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY)); 190 PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY)); 191 } 192 } 193 } else { 194 if (ts->rhsjacobianp) PetscCall(TSComputeRHSJacobianP(ts, t, U, ts->Jacprhs)); 195 if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */ 196 PetscCall(MatScale(Amat, -1)); 197 } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */ 198 MatStructure axpy = DIFFERENT_NONZERO_PATTERN; 199 if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */ 200 PetscCall(MatZeroEntries(Amat)); 201 } 202 PetscCall(MatAXPY(Amat, -1, ts->Jacprhs, axpy)); 203 } 204 } 205 PetscCall(PetscLogEventEnd(TS_JacobianPEval, ts, U, Amat, 0)); 206 PetscFunctionReturn(0); 207 } 208 209 /*@C 210 TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 211 212 Logically Collective on TS 213 214 Input Parameters: 215 + ts - the TS context obtained from TSCreate() 216 . numcost - number of gradients to be computed, this is the number of cost functions 217 . costintegral - vector that stores the integral values 218 . rf - routine for evaluating the integrand function 219 . drduf - function that computes the gradients of the r's with respect to u 220 . drdpf - function that computes the gradients of the r's with respect to p, can be NULL if parametric sensitivity is not desired (mu=NULL) 221 . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 222 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 223 224 Calling sequence of rf: 225 $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx); 226 227 Calling sequence of drduf: 228 $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx); 229 230 Calling sequence of drdpf: 231 $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx); 232 233 Level: deprecated 234 235 Notes: 236 For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 237 238 .seealso: `TSSetRHSJacobianP()`, `TSGetCostGradients()`, `TSSetCostGradients()` 239 @*/ 240 PetscErrorCode TSSetCostIntegrand(TS ts, PetscInt numcost, Vec costintegral, PetscErrorCode (*rf)(TS, PetscReal, Vec, Vec, void *), PetscErrorCode (*drduf)(TS, PetscReal, Vec, Vec *, void *), PetscErrorCode (*drdpf)(TS, PetscReal, Vec, Vec *, void *), PetscBool fwd, void *ctx) 241 { 242 PetscFunctionBegin; 243 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 244 if (costintegral) PetscValidHeaderSpecific(costintegral, VEC_CLASSID, 3); 245 PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostGradients() or TSForwardSetIntegralGradients()"); 246 if (!ts->numcost) ts->numcost = numcost; 247 248 if (costintegral) { 249 PetscCall(PetscObjectReference((PetscObject)costintegral)); 250 PetscCall(VecDestroy(&ts->vec_costintegral)); 251 ts->vec_costintegral = costintegral; 252 } else { 253 if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 254 PetscCall(VecCreateSeq(PETSC_COMM_SELF, numcost, &ts->vec_costintegral)); 255 } else { 256 PetscCall(VecSet(ts->vec_costintegral, 0.0)); 257 } 258 } 259 if (!ts->vec_costintegrand) { 260 PetscCall(VecDuplicate(ts->vec_costintegral, &ts->vec_costintegrand)); 261 } else { 262 PetscCall(VecSet(ts->vec_costintegrand, 0.0)); 263 } 264 ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 265 ts->costintegrand = rf; 266 ts->costintegrandctx = ctx; 267 ts->drdufunction = drduf; 268 ts->drdpfunction = drdpf; 269 PetscFunctionReturn(0); 270 } 271 272 /*@C 273 TSGetCostIntegral - Returns the values of the integral term in the cost functions. 274 It is valid to call the routine after a backward run. 275 276 Not Collective 277 278 Input Parameter: 279 . ts - the TS context obtained from TSCreate() 280 281 Output Parameter: 282 . v - the vector containing the integrals for each cost function 283 284 Level: intermediate 285 286 .seealso: `TSSetCostIntegrand()` 287 288 @*/ 289 PetscErrorCode TSGetCostIntegral(TS ts, Vec *v) 290 { 291 TS quadts; 292 293 PetscFunctionBegin; 294 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 295 PetscValidPointer(v, 2); 296 PetscCall(TSGetQuadratureTS(ts, NULL, &quadts)); 297 *v = quadts->vec_sol; 298 PetscFunctionReturn(0); 299 } 300 301 /*@C 302 TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 303 304 Input Parameters: 305 + ts - the TS context 306 . t - current time 307 - U - state vector, i.e. current solution 308 309 Output Parameter: 310 . Q - vector of size numcost to hold the outputs 311 312 Notes: 313 Most users should not need to explicitly call this routine, as it 314 is used internally within the sensitivity analysis context. 315 316 Level: deprecated 317 318 .seealso: `TSSetCostIntegrand()` 319 @*/ 320 PetscErrorCode TSComputeCostIntegrand(TS ts, PetscReal t, Vec U, Vec Q) 321 { 322 PetscFunctionBegin; 323 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 324 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 325 PetscValidHeaderSpecific(Q, VEC_CLASSID, 4); 326 327 PetscCall(PetscLogEventBegin(TS_FunctionEval, ts, U, Q, 0)); 328 if (ts->costintegrand) PetscCallBack("TS callback integrand in the cost function", (*ts->costintegrand)(ts, t, U, Q, ts->costintegrandctx)); 329 else PetscCall(VecZeroEntries(Q)); 330 PetscCall(PetscLogEventEnd(TS_FunctionEval, ts, U, Q, 0)); 331 PetscFunctionReturn(0); 332 } 333 334 /*@C 335 TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 336 337 Level: deprecated 338 339 @*/ 340 PetscErrorCode TSComputeDRDUFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) 341 { 342 PetscFunctionBegin; 343 if (!DRDU) PetscFunctionReturn(0); 344 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 345 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 346 347 PetscCallBack("TS callback DRDU for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 348 PetscFunctionReturn(0); 349 } 350 351 /*@C 352 TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 353 354 Level: deprecated 355 356 @*/ 357 PetscErrorCode TSComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) 358 { 359 PetscFunctionBegin; 360 if (!DRDP) PetscFunctionReturn(0); 361 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 362 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 363 364 PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 365 PetscFunctionReturn(0); 366 } 367 368 /*@C 369 TSSetIHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of F (IFunction) w.r.t. the state variable. 370 371 Logically Collective on TS 372 373 Input Parameters: 374 + ts - TS context obtained from TSCreate() 375 . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 376 . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 377 . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 378 . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 379 . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 380 . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 381 . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 382 - hessianproductfunc4 - vector-Hessian-vector product function for F_PP 383 384 Calling sequence of ihessianproductfunc: 385 $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 386 + t - current timestep 387 . U - input vector (current ODE solution) 388 . Vl - an array of input vectors to be left-multiplied with the Hessian 389 . Vr - input vector to be right-multiplied with the Hessian 390 . VHV - an array of output vectors for vector-Hessian-vector product 391 - ctx - [optional] user-defined function context 392 393 Level: intermediate 394 395 Notes: 396 The first Hessian function and the working array are required. 397 As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 398 $ Vl_n^T*F_UP*Vr 399 where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian F_UP is of size N x N x M. 400 Each entry of F_UP corresponds to the derivative 401 $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}. 402 The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with the j-th entry being 403 $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]} 404 If the cost function is a scalar, there will be only one vector in Vl and VHV. 405 406 .seealso: 407 @*/ 408 PetscErrorCode TSSetIHessianProduct(TS ts, Vec *ihp1, PetscErrorCode (*ihessianproductfunc1)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp2, PetscErrorCode (*ihessianproductfunc2)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp3, PetscErrorCode (*ihessianproductfunc3)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp4, PetscErrorCode (*ihessianproductfunc4)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *ctx) 409 { 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 412 PetscValidPointer(ihp1, 2); 413 414 ts->ihessianproductctx = ctx; 415 if (ihp1) ts->vecs_fuu = ihp1; 416 if (ihp2) ts->vecs_fup = ihp2; 417 if (ihp3) ts->vecs_fpu = ihp3; 418 if (ihp4) ts->vecs_fpp = ihp4; 419 ts->ihessianproduct_fuu = ihessianproductfunc1; 420 ts->ihessianproduct_fup = ihessianproductfunc2; 421 ts->ihessianproduct_fpu = ihessianproductfunc3; 422 ts->ihessianproduct_fpp = ihessianproductfunc4; 423 PetscFunctionReturn(0); 424 } 425 426 /*@C 427 TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu. 428 429 Collective on TS 430 431 Input Parameters: 432 . ts - The TS context obtained from TSCreate() 433 434 Notes: 435 TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation, 436 so most users would not generally call this routine themselves. 437 438 Level: developer 439 440 .seealso: `TSSetIHessianProduct()` 441 @*/ 442 PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 443 { 444 PetscFunctionBegin; 445 if (!VHV) PetscFunctionReturn(0); 446 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 447 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 448 449 if (ts->ihessianproduct_fuu) PetscCallBack("TS callback IHessianProduct 1 for sensitivity analysis", (*ts->ihessianproduct_fuu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 450 451 /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 452 if (ts->rhshessianproduct_guu) { 453 PetscInt nadj; 454 PetscCall(TSComputeRHSHessianProductFunctionUU(ts, t, U, Vl, Vr, VHV)); 455 for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 456 } 457 PetscFunctionReturn(0); 458 } 459 460 /*@C 461 TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup. 462 463 Collective on TS 464 465 Input Parameters: 466 . ts - The TS context obtained from TSCreate() 467 468 Notes: 469 TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation, 470 so most users would not generally call this routine themselves. 471 472 Level: developer 473 474 .seealso: `TSSetIHessianProduct()` 475 @*/ 476 PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 477 { 478 PetscFunctionBegin; 479 if (!VHV) PetscFunctionReturn(0); 480 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 481 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 482 483 if (ts->ihessianproduct_fup) PetscCallBack("TS callback IHessianProduct 2 for sensitivity analysis", (*ts->ihessianproduct_fup)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 484 485 /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 486 if (ts->rhshessianproduct_gup) { 487 PetscInt nadj; 488 PetscCall(TSComputeRHSHessianProductFunctionUP(ts, t, U, Vl, Vr, VHV)); 489 for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 490 } 491 PetscFunctionReturn(0); 492 } 493 494 /*@C 495 TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu. 496 497 Collective on TS 498 499 Input Parameters: 500 . ts - The TS context obtained from TSCreate() 501 502 Notes: 503 TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation, 504 so most users would not generally call this routine themselves. 505 506 Level: developer 507 508 .seealso: `TSSetIHessianProduct()` 509 @*/ 510 PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 511 { 512 PetscFunctionBegin; 513 if (!VHV) PetscFunctionReturn(0); 514 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 515 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 516 517 if (ts->ihessianproduct_fpu) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 518 519 /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 520 if (ts->rhshessianproduct_gpu) { 521 PetscInt nadj; 522 PetscCall(TSComputeRHSHessianProductFunctionPU(ts, t, U, Vl, Vr, VHV)); 523 for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 524 } 525 PetscFunctionReturn(0); 526 } 527 528 /*@C 529 TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp. 530 531 Collective on TS 532 533 Input Parameters: 534 . ts - The TS context obtained from TSCreate() 535 536 Notes: 537 TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation, 538 so most users would not generally call this routine themselves. 539 540 Level: developer 541 542 .seealso: `TSSetIHessianProduct()` 543 @*/ 544 PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 545 { 546 PetscFunctionBegin; 547 if (!VHV) PetscFunctionReturn(0); 548 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 549 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 550 551 if (ts->ihessianproduct_fpp) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpp)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 552 553 /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 554 if (ts->rhshessianproduct_gpp) { 555 PetscInt nadj; 556 PetscCall(TSComputeRHSHessianProductFunctionPP(ts, t, U, Vl, Vr, VHV)); 557 for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 558 } 559 PetscFunctionReturn(0); 560 } 561 562 /*@C 563 TSSetRHSHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of G (RHSFunction) w.r.t. the state variable. 564 565 Logically Collective on TS 566 567 Input Parameters: 568 + ts - TS context obtained from TSCreate() 569 . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 570 . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 571 . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 572 . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 573 . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 574 . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 575 . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 576 - hessianproductfunc4 - vector-Hessian-vector product function for G_PP 577 578 Calling sequence of ihessianproductfunc: 579 $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 580 + t - current timestep 581 . U - input vector (current ODE solution) 582 . Vl - an array of input vectors to be left-multiplied with the Hessian 583 . Vr - input vector to be right-multiplied with the Hessian 584 . VHV - an array of output vectors for vector-Hessian-vector product 585 - ctx - [optional] user-defined function context 586 587 Level: intermediate 588 589 Notes: 590 The first Hessian function and the working array are required. 591 As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 592 $ Vl_n^T*G_UP*Vr 593 where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian G_UP is of size N x N x M. 594 Each entry of G_UP corresponds to the derivative 595 $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}. 596 The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with j-th entry being 597 $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]} 598 If the cost function is a scalar, there will be only one vector in Vl and VHV. 599 600 .seealso: 601 @*/ 602 PetscErrorCode TSSetRHSHessianProduct(TS ts, Vec *rhshp1, PetscErrorCode (*rhshessianproductfunc1)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp2, PetscErrorCode (*rhshessianproductfunc2)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp3, PetscErrorCode (*rhshessianproductfunc3)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp4, PetscErrorCode (*rhshessianproductfunc4)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *ctx) 603 { 604 PetscFunctionBegin; 605 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 606 PetscValidPointer(rhshp1, 2); 607 608 ts->rhshessianproductctx = ctx; 609 if (rhshp1) ts->vecs_guu = rhshp1; 610 if (rhshp2) ts->vecs_gup = rhshp2; 611 if (rhshp3) ts->vecs_gpu = rhshp3; 612 if (rhshp4) ts->vecs_gpp = rhshp4; 613 ts->rhshessianproduct_guu = rhshessianproductfunc1; 614 ts->rhshessianproduct_gup = rhshessianproductfunc2; 615 ts->rhshessianproduct_gpu = rhshessianproductfunc3; 616 ts->rhshessianproduct_gpp = rhshessianproductfunc4; 617 PetscFunctionReturn(0); 618 } 619 620 /*@C 621 TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu. 622 623 Collective on TS 624 625 Input Parameters: 626 . ts - The TS context obtained from TSCreate() 627 628 Notes: 629 TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation, 630 so most users would not generally call this routine themselves. 631 632 Level: developer 633 634 .seealso: `TSSetRHSHessianProduct()` 635 @*/ 636 PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 637 { 638 PetscFunctionBegin; 639 if (!VHV) PetscFunctionReturn(0); 640 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 641 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 642 643 PetscCallBack("TS callback RHSHessianProduct 1 for sensitivity analysis", (*ts->rhshessianproduct_guu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 644 PetscFunctionReturn(0); 645 } 646 647 /*@C 648 TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup. 649 650 Collective on TS 651 652 Input Parameters: 653 . ts - The TS context obtained from TSCreate() 654 655 Notes: 656 TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation, 657 so most users would not generally call this routine themselves. 658 659 Level: developer 660 661 .seealso: `TSSetRHSHessianProduct()` 662 @*/ 663 PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 664 { 665 PetscFunctionBegin; 666 if (!VHV) PetscFunctionReturn(0); 667 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 668 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 669 670 PetscCallBack("TS callback RHSHessianProduct 2 for sensitivity analysis", (*ts->rhshessianproduct_gup)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 671 PetscFunctionReturn(0); 672 } 673 674 /*@C 675 TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu. 676 677 Collective on TS 678 679 Input Parameters: 680 . ts - The TS context obtained from TSCreate() 681 682 Notes: 683 TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation, 684 so most users would not generally call this routine themselves. 685 686 Level: developer 687 688 .seealso: `TSSetRHSHessianProduct()` 689 @*/ 690 PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 691 { 692 PetscFunctionBegin; 693 if (!VHV) PetscFunctionReturn(0); 694 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 695 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 696 697 PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 698 PetscFunctionReturn(0); 699 } 700 701 /*@C 702 TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp. 703 704 Collective on TS 705 706 Input Parameters: 707 . ts - The TS context obtained from TSCreate() 708 709 Notes: 710 TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation, 711 so most users would not generally call this routine themselves. 712 713 Level: developer 714 715 .seealso: `TSSetRHSHessianProduct()` 716 @*/ 717 PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 718 { 719 PetscFunctionBegin; 720 if (!VHV) PetscFunctionReturn(0); 721 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 722 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 723 724 PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpp)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 725 PetscFunctionReturn(0); 726 } 727 728 /* --------------------------- Adjoint sensitivity ---------------------------*/ 729 730 /*@ 731 TSSetCostGradients - Sets the initial value of the gradients of the cost function w.r.t. initial values and w.r.t. the problem parameters 732 for use by the TSAdjoint routines. 733 734 Logically Collective on TS 735 736 Input Parameters: 737 + ts - the TS context obtained from TSCreate() 738 . numcost - number of gradients to be computed, this is the number of cost functions 739 . lambda - gradients with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector 740 - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 741 742 Level: beginner 743 744 Notes: 745 the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 746 747 After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 748 749 .seealso `TSGetCostGradients()` 750 @*/ 751 PetscErrorCode TSSetCostGradients(TS ts, PetscInt numcost, Vec *lambda, Vec *mu) 752 { 753 PetscFunctionBegin; 754 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 755 PetscValidPointer(lambda, 3); 756 ts->vecs_sensi = lambda; 757 ts->vecs_sensip = mu; 758 PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand"); 759 ts->numcost = numcost; 760 PetscFunctionReturn(0); 761 } 762 763 /*@ 764 TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 765 766 Not Collective, but Vec returned is parallel if TS is parallel 767 768 Input Parameter: 769 . ts - the TS context obtained from TSCreate() 770 771 Output Parameters: 772 + numcost - size of returned arrays 773 . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 774 - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 775 776 Level: intermediate 777 778 .seealso: `TSSetCostGradients()` 779 @*/ 780 PetscErrorCode TSGetCostGradients(TS ts, PetscInt *numcost, Vec **lambda, Vec **mu) 781 { 782 PetscFunctionBegin; 783 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 784 if (numcost) *numcost = ts->numcost; 785 if (lambda) *lambda = ts->vecs_sensi; 786 if (mu) *mu = ts->vecs_sensip; 787 PetscFunctionReturn(0); 788 } 789 790 /*@ 791 TSSetCostHessianProducts - Sets the initial value of the Hessian-vector products of the cost function w.r.t. initial values and w.r.t. the problem parameters 792 for use by the TSAdjoint routines. 793 794 Logically Collective on TS 795 796 Input Parameters: 797 + ts - the TS context obtained from TSCreate() 798 . numcost - number of cost functions 799 . lambda2 - Hessian-vector product with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector 800 . mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 801 - dir - the direction vector that are multiplied with the Hessian of the cost functions 802 803 Level: beginner 804 805 Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system 806 807 For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve(). 808 809 After TSAdjointSolve() is called, the lamba2 and the mu2 will contain the computed second-order adjoint sensitivities, and can be used to produce Hessian-vector product (not the full Hessian matrix). Users must provide a direction vector; it is usually generated by an optimization solver. 810 811 Passing NULL for lambda2 disables the second-order calculation. 812 .seealso: `TSAdjointSetForward()` 813 @*/ 814 PetscErrorCode TSSetCostHessianProducts(TS ts, PetscInt numcost, Vec *lambda2, Vec *mu2, Vec dir) 815 { 816 PetscFunctionBegin; 817 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 818 PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand"); 819 ts->numcost = numcost; 820 ts->vecs_sensi2 = lambda2; 821 ts->vecs_sensi2p = mu2; 822 ts->vec_dir = dir; 823 PetscFunctionReturn(0); 824 } 825 826 /*@ 827 TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve() 828 829 Not Collective, but Vec returned is parallel if TS is parallel 830 831 Input Parameter: 832 . ts - the TS context obtained from TSCreate() 833 834 Output Parameters: 835 + numcost - number of cost functions 836 . lambda2 - Hessian-vector product with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector 837 . mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 838 - dir - the direction vector that are multiplied with the Hessian of the cost functions 839 840 Level: intermediate 841 842 .seealso: `TSSetCostHessianProducts()` 843 @*/ 844 PetscErrorCode TSGetCostHessianProducts(TS ts, PetscInt *numcost, Vec **lambda2, Vec **mu2, Vec *dir) 845 { 846 PetscFunctionBegin; 847 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 848 if (numcost) *numcost = ts->numcost; 849 if (lambda2) *lambda2 = ts->vecs_sensi2; 850 if (mu2) *mu2 = ts->vecs_sensi2p; 851 if (dir) *dir = ts->vec_dir; 852 PetscFunctionReturn(0); 853 } 854 855 /*@ 856 TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities 857 858 Logically Collective on TS 859 860 Input Parameters: 861 + ts - the TS context obtained from TSCreate() 862 - didp - the derivative of initial values w.r.t. parameters 863 864 Level: intermediate 865 866 Notes: When computing sensitivies w.r.t. initial condition, set didp to NULL so that the solver will take it as an identity matrix mathematically. TSAdjoint does not reset the tangent linear solver automatically, TSAdjointResetForward() should be called to reset the tangent linear solver. 867 868 .seealso: `TSSetCostHessianProducts()`, `TSAdjointResetForward()` 869 @*/ 870 PetscErrorCode TSAdjointSetForward(TS ts, Mat didp) 871 { 872 Mat A; 873 Vec sp; 874 PetscScalar *xarr; 875 PetscInt lsize; 876 877 PetscFunctionBegin; 878 ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 879 PetscCheck(ts->vecs_sensi2, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetCostHessianProducts() first"); 880 PetscCheck(ts->vec_dir, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 881 /* create a single-column dense matrix */ 882 PetscCall(VecGetLocalSize(ts->vec_sol, &lsize)); 883 PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ts), lsize, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, &A)); 884 885 PetscCall(VecDuplicate(ts->vec_sol, &sp)); 886 PetscCall(MatDenseGetColumn(A, 0, &xarr)); 887 PetscCall(VecPlaceArray(sp, xarr)); 888 if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */ 889 if (didp) { 890 PetscCall(MatMult(didp, ts->vec_dir, sp)); 891 PetscCall(VecScale(sp, 2.)); 892 } else { 893 PetscCall(VecZeroEntries(sp)); 894 } 895 } else { /* tangent linear variable initialized as dir */ 896 PetscCall(VecCopy(ts->vec_dir, sp)); 897 } 898 PetscCall(VecResetArray(sp)); 899 PetscCall(MatDenseRestoreColumn(A, &xarr)); 900 PetscCall(VecDestroy(&sp)); 901 902 PetscCall(TSForwardSetInitialSensitivities(ts, A)); /* if didp is NULL, identity matrix is assumed */ 903 904 PetscCall(MatDestroy(&A)); 905 PetscFunctionReturn(0); 906 } 907 908 /*@ 909 TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context 910 911 Logically Collective on TS 912 913 Input Parameters: 914 . ts - the TS context obtained from TSCreate() 915 916 Level: intermediate 917 918 .seealso: `TSAdjointSetForward()` 919 @*/ 920 PetscErrorCode TSAdjointResetForward(TS ts) 921 { 922 PetscFunctionBegin; 923 ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */ 924 PetscCall(TSForwardReset(ts)); 925 PetscFunctionReturn(0); 926 } 927 928 /*@ 929 TSAdjointSetUp - Sets up the internal data structures for the later use 930 of an adjoint solver 931 932 Collective on TS 933 934 Input Parameter: 935 . ts - the TS context obtained from TSCreate() 936 937 Level: advanced 938 939 .seealso: `TSCreate()`, `TSAdjointStep()`, `TSSetCostGradients()` 940 @*/ 941 PetscErrorCode TSAdjointSetUp(TS ts) 942 { 943 TSTrajectory tj; 944 PetscBool match; 945 946 PetscFunctionBegin; 947 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 948 if (ts->adjointsetupcalled) PetscFunctionReturn(0); 949 PetscCheck(ts->vecs_sensi, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetCostGradients() first"); 950 PetscCheck(!ts->vecs_sensip || ts->Jacp || ts->Jacprhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 951 PetscCall(TSGetTrajectory(ts, &tj)); 952 PetscCall(PetscObjectTypeCompare((PetscObject)tj, TSTRAJECTORYBASIC, &match)); 953 if (match) { 954 PetscBool solution_only; 955 PetscCall(TSTrajectoryGetSolutionOnly(tj, &solution_only)); 956 PetscCheck(!solution_only, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "TSAdjoint cannot use the solution-only mode when choosing the Basic TSTrajectory type. Turn it off with -ts_trajectory_solution_only 0"); 957 } 958 PetscCall(TSTrajectorySetUseHistory(tj, PETSC_FALSE)); /* not use TSHistory */ 959 960 if (ts->quadraturets) { /* if there is integral in the cost function */ 961 PetscCall(VecDuplicate(ts->vecs_sensi[0], &ts->vec_drdu_col)); 962 if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &ts->vec_drdp_col)); 963 } 964 965 PetscTryTypeMethod(ts, adjointsetup); 966 ts->adjointsetupcalled = PETSC_TRUE; 967 PetscFunctionReturn(0); 968 } 969 970 /*@ 971 TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 972 973 Collective on TS 974 975 Input Parameter: 976 . ts - the TS context obtained from TSCreate() 977 978 Level: beginner 979 980 .seealso: `TSCreate()`, `TSAdjointSetUp()`, `TSADestroy()` 981 @*/ 982 PetscErrorCode TSAdjointReset(TS ts) 983 { 984 PetscFunctionBegin; 985 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 986 PetscTryTypeMethod(ts, adjointreset); 987 if (ts->quadraturets) { /* if there is integral in the cost function */ 988 PetscCall(VecDestroy(&ts->vec_drdu_col)); 989 if (ts->vecs_sensip) PetscCall(VecDestroy(&ts->vec_drdp_col)); 990 } 991 ts->vecs_sensi = NULL; 992 ts->vecs_sensip = NULL; 993 ts->vecs_sensi2 = NULL; 994 ts->vecs_sensi2p = NULL; 995 ts->vec_dir = NULL; 996 ts->adjointsetupcalled = PETSC_FALSE; 997 PetscFunctionReturn(0); 998 } 999 1000 /*@ 1001 TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1002 1003 Logically Collective on TS 1004 1005 Input Parameters: 1006 + ts - the TS context obtained from TSCreate() 1007 - steps - number of steps to use 1008 1009 Level: intermediate 1010 1011 Notes: 1012 Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 1013 so as to integrate back to less than the original timestep 1014 1015 .seealso: `TSSetExactFinalTime()` 1016 @*/ 1017 PetscErrorCode TSAdjointSetSteps(TS ts, PetscInt steps) 1018 { 1019 PetscFunctionBegin; 1020 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1021 PetscValidLogicalCollectiveInt(ts, steps, 2); 1022 PetscCheck(steps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back a negative number of steps"); 1023 PetscCheck(steps <= ts->steps, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back more than the total number of forward steps"); 1024 ts->adjoint_max_steps = steps; 1025 PetscFunctionReturn(0); 1026 } 1027 1028 /*@C 1029 TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1030 1031 Level: deprecated 1032 1033 @*/ 1034 PetscErrorCode TSAdjointSetRHSJacobian(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) 1035 { 1036 PetscFunctionBegin; 1037 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1038 PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 1039 1040 ts->rhsjacobianp = func; 1041 ts->rhsjacobianpctx = ctx; 1042 if (Amat) { 1043 PetscCall(PetscObjectReference((PetscObject)Amat)); 1044 PetscCall(MatDestroy(&ts->Jacp)); 1045 ts->Jacp = Amat; 1046 } 1047 PetscFunctionReturn(0); 1048 } 1049 1050 /*@C 1051 TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1052 1053 Level: deprecated 1054 1055 @*/ 1056 PetscErrorCode TSAdjointComputeRHSJacobian(TS ts, PetscReal t, Vec U, Mat Amat) 1057 { 1058 PetscFunctionBegin; 1059 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1060 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1061 PetscValidPointer(Amat, 4); 1062 1063 PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 1064 PetscFunctionReturn(0); 1065 } 1066 1067 /*@ 1068 TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 1069 1070 Level: deprecated 1071 1072 @*/ 1073 PetscErrorCode TSAdjointComputeDRDYFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) 1074 { 1075 PetscFunctionBegin; 1076 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1077 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1078 1079 PetscCallBack("TS callback DRDY for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 1080 PetscFunctionReturn(0); 1081 } 1082 1083 /*@ 1084 TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 1085 1086 Level: deprecated 1087 1088 @*/ 1089 PetscErrorCode TSAdjointComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) 1090 { 1091 PetscFunctionBegin; 1092 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1093 PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1094 1095 PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 1096 PetscFunctionReturn(0); 1097 } 1098 1099 /*@C 1100 TSAdjointMonitorSensi - monitors the first lambda sensitivity 1101 1102 Level: intermediate 1103 1104 .seealso: `TSAdjointMonitorSet()` 1105 @*/ 1106 PetscErrorCode TSAdjointMonitorSensi(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) 1107 { 1108 PetscViewer viewer = vf->viewer; 1109 1110 PetscFunctionBegin; 1111 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 1112 PetscCall(PetscViewerPushFormat(viewer, vf->format)); 1113 PetscCall(VecView(lambda[0], viewer)); 1114 PetscCall(PetscViewerPopFormat(viewer)); 1115 PetscFunctionReturn(0); 1116 } 1117 1118 /*@C 1119 TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1120 1121 Collective on TS 1122 1123 Input Parameters: 1124 + ts - TS object you wish to monitor 1125 . name - the monitor type one is seeking 1126 . help - message indicating what monitoring is done 1127 . manual - manual page for the monitor 1128 . monitor - the monitor function 1129 - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects 1130 1131 Level: developer 1132 1133 .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 1134 `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 1135 `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1136 `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1137 `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1138 `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1139 `PetscOptionsFList()`, `PetscOptionsEList()` 1140 @*/ 1141 PetscErrorCode TSAdjointMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *)) 1142 { 1143 PetscViewer viewer; 1144 PetscViewerFormat format; 1145 PetscBool flg; 1146 1147 PetscFunctionBegin; 1148 PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg)); 1149 if (flg) { 1150 PetscViewerAndFormat *vf; 1151 PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 1152 PetscCall(PetscObjectDereference((PetscObject)viewer)); 1153 if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 1154 PetscCall(TSAdjointMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 1155 } 1156 PetscFunctionReturn(0); 1157 } 1158 1159 /*@C 1160 TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1161 timestep to display the iteration's progress. 1162 1163 Logically Collective on TS 1164 1165 Input Parameters: 1166 + ts - the TS context obtained from TSCreate() 1167 . adjointmonitor - monitoring routine 1168 . adjointmctx - [optional] user-defined context for private data for the 1169 monitor routine (use NULL if no context is desired) 1170 - adjointmonitordestroy - [optional] routine that frees monitor context 1171 (may be NULL) 1172 1173 Calling sequence of monitor: 1174 $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1175 1176 + ts - the TS context 1177 . steps - iteration number (after the final time step the monitor routine is called with a step of -1, this is at the final time which may have 1178 been interpolated to) 1179 . time - current time 1180 . u - current iterate 1181 . numcost - number of cost functionos 1182 . lambda - sensitivities to initial conditions 1183 . mu - sensitivities to parameters 1184 - adjointmctx - [optional] adjoint monitoring context 1185 1186 Notes: 1187 This routine adds an additional monitor to the list of monitors that 1188 already has been loaded. 1189 1190 Fortran Notes: 1191 Only a single monitor function can be set for each TS object 1192 1193 Level: intermediate 1194 1195 .seealso: `TSAdjointMonitorCancel()` 1196 @*/ 1197 PetscErrorCode TSAdjointMonitorSet(TS ts, PetscErrorCode (*adjointmonitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *adjointmctx, PetscErrorCode (*adjointmdestroy)(void **)) 1198 { 1199 PetscInt i; 1200 PetscBool identical; 1201 1202 PetscFunctionBegin; 1203 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1204 for (i = 0; i < ts->numbermonitors; i++) { 1205 PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))adjointmonitor, adjointmctx, adjointmdestroy, (PetscErrorCode(*)(void))ts->adjointmonitor[i], ts->adjointmonitorcontext[i], ts->adjointmonitordestroy[i], &identical)); 1206 if (identical) PetscFunctionReturn(0); 1207 } 1208 PetscCheck(ts->numberadjointmonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many adjoint monitors set"); 1209 ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1210 ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1211 ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void *)adjointmctx; 1212 PetscFunctionReturn(0); 1213 } 1214 1215 /*@C 1216 TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1217 1218 Logically Collective on TS 1219 1220 Input Parameters: 1221 . ts - the TS context obtained from TSCreate() 1222 1223 Notes: 1224 There is no way to remove a single, specific monitor. 1225 1226 Level: intermediate 1227 1228 .seealso: `TSAdjointMonitorSet()` 1229 @*/ 1230 PetscErrorCode TSAdjointMonitorCancel(TS ts) 1231 { 1232 PetscInt i; 1233 1234 PetscFunctionBegin; 1235 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1236 for (i = 0; i < ts->numberadjointmonitors; i++) { 1237 if (ts->adjointmonitordestroy[i]) PetscCall((*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i])); 1238 } 1239 ts->numberadjointmonitors = 0; 1240 PetscFunctionReturn(0); 1241 } 1242 1243 /*@C 1244 TSAdjointMonitorDefault - the default monitor of adjoint computations 1245 1246 Level: intermediate 1247 1248 .seealso: `TSAdjointMonitorSet()` 1249 @*/ 1250 PetscErrorCode TSAdjointMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) 1251 { 1252 PetscViewer viewer = vf->viewer; 1253 1254 PetscFunctionBegin; 1255 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 1256 PetscCall(PetscViewerPushFormat(viewer, vf->format)); 1257 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 1258 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 1259 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 1260 PetscCall(PetscViewerPopFormat(viewer)); 1261 PetscFunctionReturn(0); 1262 } 1263 1264 /*@C 1265 TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1266 VecView() for the sensitivities to initial states at each timestep 1267 1268 Collective on TS 1269 1270 Input Parameters: 1271 + ts - the TS context 1272 . step - current time-step 1273 . ptime - current time 1274 . u - current state 1275 . numcost - number of cost functions 1276 . lambda - sensitivities to initial conditions 1277 . mu - sensitivities to parameters 1278 - dummy - either a viewer or NULL 1279 1280 Level: intermediate 1281 1282 .seealso: `TSAdjointMonitorSet()`, `TSAdjointMonitorDefault()`, `VecView()` 1283 @*/ 1284 PetscErrorCode TSAdjointMonitorDrawSensi(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *dummy) 1285 { 1286 TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1287 PetscDraw draw; 1288 PetscReal xl, yl, xr, yr, h; 1289 char time[32]; 1290 1291 PetscFunctionBegin; 1292 if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1293 1294 PetscCall(VecView(lambda[0], ictx->viewer)); 1295 PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 1296 PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 1297 PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 1298 h = yl + .95 * (yr - yl); 1299 PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 1300 PetscCall(PetscDrawFlush(draw)); 1301 PetscFunctionReturn(0); 1302 } 1303 1304 /* 1305 TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1306 1307 Collective on TSAdjoint 1308 1309 Input Parameter: 1310 . ts - the TS context 1311 1312 Options Database Keys: 1313 + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1314 . -ts_adjoint_monitor - print information at each adjoint time step 1315 - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1316 1317 Level: developer 1318 1319 Notes: 1320 This is not normally called directly by users 1321 1322 .seealso: `TSSetSaveTrajectory()`, `TSTrajectorySetUp()` 1323 */ 1324 PetscErrorCode TSAdjointSetFromOptions(TS ts, PetscOptionItems *PetscOptionsObject) 1325 { 1326 PetscBool tflg, opt; 1327 1328 PetscFunctionBegin; 1329 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1330 PetscOptionsHeadBegin(PetscOptionsObject, "TS Adjoint options"); 1331 tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 1332 PetscCall(PetscOptionsBool("-ts_adjoint_solve", "Solve the adjoint problem immediately after solving the forward problem", "", tflg, &tflg, &opt)); 1333 if (opt) { 1334 PetscCall(TSSetSaveTrajectory(ts)); 1335 ts->adjoint_solve = tflg; 1336 } 1337 PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor", "Monitor adjoint timestep size", "TSAdjointMonitorDefault", TSAdjointMonitorDefault, NULL)); 1338 PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor_sensi", "Monitor sensitivity in the adjoint computation", "TSAdjointMonitorSensi", TSAdjointMonitorSensi, NULL)); 1339 opt = PETSC_FALSE; 1340 PetscCall(PetscOptionsName("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", &opt)); 1341 if (opt) { 1342 TSMonitorDrawCtx ctx; 1343 PetscInt howoften = 1; 1344 1345 PetscCall(PetscOptionsInt("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", howoften, &howoften, NULL)); 1346 PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx)); 1347 PetscCall(TSAdjointMonitorSet(ts, TSAdjointMonitorDrawSensi, ctx, (PetscErrorCode(*)(void **))TSMonitorDrawCtxDestroy)); 1348 } 1349 PetscFunctionReturn(0); 1350 } 1351 1352 /*@ 1353 TSAdjointStep - Steps one time step backward in the adjoint run 1354 1355 Collective on TS 1356 1357 Input Parameter: 1358 . ts - the TS context obtained from TSCreate() 1359 1360 Level: intermediate 1361 1362 .seealso: `TSAdjointSetUp()`, `TSAdjointSolve()` 1363 @*/ 1364 PetscErrorCode TSAdjointStep(TS ts) 1365 { 1366 DM dm; 1367 1368 PetscFunctionBegin; 1369 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1370 PetscCall(TSGetDM(ts, &dm)); 1371 PetscCall(TSAdjointSetUp(ts)); 1372 ts->steps--; /* must decrease the step index before the adjoint step is taken. */ 1373 1374 ts->reason = TS_CONVERGED_ITERATING; 1375 ts->ptime_prev = ts->ptime; 1376 PetscCall(PetscLogEventBegin(TS_AdjointStep, ts, 0, 0, 0)); 1377 PetscUseTypeMethod(ts, adjointstep); 1378 PetscCall(PetscLogEventEnd(TS_AdjointStep, ts, 0, 0, 0)); 1379 ts->adjoint_steps++; 1380 1381 if (ts->reason < 0) { 1382 PetscCheck(!ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSAdjointStep has failed due to %s", TSConvergedReasons[ts->reason]); 1383 } else if (!ts->reason) { 1384 if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1385 } 1386 PetscFunctionReturn(0); 1387 } 1388 1389 /*@ 1390 TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1391 1392 Collective on TS 1393 1394 Input Parameter: 1395 . ts - the TS context obtained from TSCreate() 1396 1397 Options Database: 1398 . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1399 1400 Level: intermediate 1401 1402 Notes: 1403 This must be called after a call to TSSolve() that solves the forward problem 1404 1405 By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1406 1407 .seealso: `TSCreate()`, `TSSetCostGradients()`, `TSSetSolution()`, `TSAdjointStep()` 1408 @*/ 1409 PetscErrorCode TSAdjointSolve(TS ts) 1410 { 1411 static PetscBool cite = PETSC_FALSE; 1412 #if defined(TSADJOINT_STAGE) 1413 PetscLogStage adjoint_stage; 1414 #endif 1415 1416 PetscFunctionBegin; 1417 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1418 PetscCall(PetscCitationsRegister("@article{Zhang2022tsadjoint,\n" 1419 " title = {{PETSc TSAdjoint: A Discrete Adjoint ODE Solver for First-Order and Second-Order Sensitivity Analysis}},\n" 1420 " author = {Zhang, Hong and Constantinescu, Emil M. and Smith, Barry F.},\n" 1421 " journal = {SIAM Journal on Scientific Computing},\n" 1422 " volume = {44},\n" 1423 " number = {1},\n" 1424 " pages = {C1-C24},\n" 1425 " doi = {10.1137/21M140078X},\n" 1426 " year = {2022}\n}\n", 1427 &cite)); 1428 #if defined(TSADJOINT_STAGE) 1429 PetscCall(PetscLogStageRegister("TSAdjoint", &adjoint_stage)); 1430 PetscCall(PetscLogStagePush(adjoint_stage)); 1431 #endif 1432 PetscCall(TSAdjointSetUp(ts)); 1433 1434 /* reset time step and iteration counters */ 1435 ts->adjoint_steps = 0; 1436 ts->ksp_its = 0; 1437 ts->snes_its = 0; 1438 ts->num_snes_failures = 0; 1439 ts->reject = 0; 1440 ts->reason = TS_CONVERGED_ITERATING; 1441 1442 if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1443 if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1444 1445 while (!ts->reason) { 1446 PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 1447 PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 1448 PetscCall(TSAdjointEventHandler(ts)); 1449 PetscCall(TSAdjointStep(ts)); 1450 if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) PetscCall(TSAdjointCostIntegral(ts)); 1451 } 1452 if (!ts->steps) { 1453 PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 1454 PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 1455 } 1456 ts->solvetime = ts->ptime; 1457 PetscCall(TSTrajectoryViewFromOptions(ts->trajectory, NULL, "-ts_trajectory_view")); 1458 PetscCall(VecViewFromOptions(ts->vecs_sensi[0], (PetscObject)ts, "-ts_adjoint_view_solution")); 1459 ts->adjoint_max_steps = 0; 1460 #if defined(TSADJOINT_STAGE) 1461 PetscCall(PetscLogStagePop()); 1462 #endif 1463 PetscFunctionReturn(0); 1464 } 1465 1466 /*@C 1467 TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1468 1469 Collective on TS 1470 1471 Input Parameters: 1472 + ts - time stepping context obtained from TSCreate() 1473 . step - step number that has just completed 1474 . ptime - model time of the state 1475 . u - state at the current model time 1476 . numcost - number of cost functions (dimension of lambda or mu) 1477 . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1478 - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1479 1480 Notes: 1481 TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1482 Users would almost never call this routine directly. 1483 1484 Level: developer 1485 1486 @*/ 1487 PetscErrorCode TSAdjointMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu) 1488 { 1489 PetscInt i, n = ts->numberadjointmonitors; 1490 1491 PetscFunctionBegin; 1492 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1493 PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 1494 PetscCall(VecLockReadPush(u)); 1495 for (i = 0; i < n; i++) PetscCall((*ts->adjointmonitor[i])(ts, step, ptime, u, numcost, lambda, mu, ts->adjointmonitorcontext[i])); 1496 PetscCall(VecLockReadPop(u)); 1497 PetscFunctionReturn(0); 1498 } 1499 1500 /*@ 1501 TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1502 1503 Collective on TS 1504 1505 Input Parameter: 1506 . ts - time stepping context 1507 1508 Level: advanced 1509 1510 Notes: 1511 This function cannot be called until TSAdjointStep() has been completed. 1512 1513 .seealso: `TSAdjointSolve()`, `TSAdjointStep` 1514 @*/ 1515 PetscErrorCode TSAdjointCostIntegral(TS ts) 1516 { 1517 PetscFunctionBegin; 1518 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1519 PetscUseTypeMethod(ts, adjointintegral); 1520 PetscFunctionReturn(0); 1521 } 1522 1523 /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1524 1525 /*@ 1526 TSForwardSetUp - Sets up the internal data structures for the later use 1527 of forward sensitivity analysis 1528 1529 Collective on TS 1530 1531 Input Parameter: 1532 . ts - the TS context obtained from TSCreate() 1533 1534 Level: advanced 1535 1536 .seealso: `TSCreate()`, `TSDestroy()`, `TSSetUp()` 1537 @*/ 1538 PetscErrorCode TSForwardSetUp(TS ts) 1539 { 1540 PetscFunctionBegin; 1541 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1542 if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1543 PetscTryTypeMethod(ts, forwardsetup); 1544 PetscCall(VecDuplicate(ts->vec_sol, &ts->vec_sensip_col)); 1545 ts->forwardsetupcalled = PETSC_TRUE; 1546 PetscFunctionReturn(0); 1547 } 1548 1549 /*@ 1550 TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 1551 1552 Collective on TS 1553 1554 Input Parameter: 1555 . ts - the TS context obtained from TSCreate() 1556 1557 Level: advanced 1558 1559 .seealso: `TSCreate()`, `TSDestroy()`, `TSForwardSetUp()` 1560 @*/ 1561 PetscErrorCode TSForwardReset(TS ts) 1562 { 1563 TS quadts = ts->quadraturets; 1564 1565 PetscFunctionBegin; 1566 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1567 PetscTryTypeMethod(ts, forwardreset); 1568 PetscCall(MatDestroy(&ts->mat_sensip)); 1569 if (quadts) PetscCall(MatDestroy(&quadts->mat_sensip)); 1570 PetscCall(VecDestroy(&ts->vec_sensip_col)); 1571 ts->forward_solve = PETSC_FALSE; 1572 ts->forwardsetupcalled = PETSC_FALSE; 1573 PetscFunctionReturn(0); 1574 } 1575 1576 /*@ 1577 TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1578 1579 Input Parameters: 1580 + ts - the TS context obtained from TSCreate() 1581 . numfwdint - number of integrals 1582 - vp - the vectors containing the gradients for each integral w.r.t. parameters 1583 1584 Level: deprecated 1585 1586 .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1587 @*/ 1588 PetscErrorCode TSForwardSetIntegralGradients(TS ts, PetscInt numfwdint, Vec *vp) 1589 { 1590 PetscFunctionBegin; 1591 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1592 PetscCheck(!ts->numcost || ts->numcost == numfwdint, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand()"); 1593 if (!ts->numcost) ts->numcost = numfwdint; 1594 1595 ts->vecs_integral_sensip = vp; 1596 PetscFunctionReturn(0); 1597 } 1598 1599 /*@ 1600 TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1601 1602 Input Parameter: 1603 . ts - the TS context obtained from TSCreate() 1604 1605 Output Parameter: 1606 . vp - the vectors containing the gradients for each integral w.r.t. parameters 1607 1608 Level: deprecated 1609 1610 .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1611 @*/ 1612 PetscErrorCode TSForwardGetIntegralGradients(TS ts, PetscInt *numfwdint, Vec **vp) 1613 { 1614 PetscFunctionBegin; 1615 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1616 PetscValidPointer(vp, 3); 1617 if (numfwdint) *numfwdint = ts->numcost; 1618 if (vp) *vp = ts->vecs_integral_sensip; 1619 PetscFunctionReturn(0); 1620 } 1621 1622 /*@ 1623 TSForwardStep - Compute the forward sensitivity for one time step. 1624 1625 Collective on TS 1626 1627 Input Parameter: 1628 . ts - time stepping context 1629 1630 Level: advanced 1631 1632 Notes: 1633 This function cannot be called until TSStep() has been completed. 1634 1635 .seealso: `TSForwardSetSensitivities()`, `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardSetUp()` 1636 @*/ 1637 PetscErrorCode TSForwardStep(TS ts) 1638 { 1639 PetscFunctionBegin; 1640 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1641 PetscCall(PetscLogEventBegin(TS_ForwardStep, ts, 0, 0, 0)); 1642 PetscUseTypeMethod(ts, forwardstep); 1643 PetscCall(PetscLogEventEnd(TS_ForwardStep, ts, 0, 0, 0)); 1644 PetscCheck(ts->reason >= 0 || !ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSFowardStep has failed due to %s", TSConvergedReasons[ts->reason]); 1645 PetscFunctionReturn(0); 1646 } 1647 1648 /*@ 1649 TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1650 1651 Logically Collective on TS 1652 1653 Input Parameters: 1654 + ts - the TS context obtained from TSCreate() 1655 . nump - number of parameters 1656 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1657 1658 Level: beginner 1659 1660 Notes: 1661 Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1662 This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1663 You must call this function before TSSolve(). 1664 The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1665 1666 .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1667 @*/ 1668 PetscErrorCode TSForwardSetSensitivities(TS ts, PetscInt nump, Mat Smat) 1669 { 1670 PetscFunctionBegin; 1671 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1672 PetscValidHeaderSpecific(Smat, MAT_CLASSID, 3); 1673 ts->forward_solve = PETSC_TRUE; 1674 if (nump == PETSC_DEFAULT) { 1675 PetscCall(MatGetSize(Smat, NULL, &ts->num_parameters)); 1676 } else ts->num_parameters = nump; 1677 PetscCall(PetscObjectReference((PetscObject)Smat)); 1678 PetscCall(MatDestroy(&ts->mat_sensip)); 1679 ts->mat_sensip = Smat; 1680 PetscFunctionReturn(0); 1681 } 1682 1683 /*@ 1684 TSForwardGetSensitivities - Returns the trajectory sensitivities 1685 1686 Not Collective, but Vec returned is parallel if TS is parallel 1687 1688 Output Parameters: 1689 + ts - the TS context obtained from TSCreate() 1690 . nump - number of parameters 1691 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1692 1693 Level: intermediate 1694 1695 .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1696 @*/ 1697 PetscErrorCode TSForwardGetSensitivities(TS ts, PetscInt *nump, Mat *Smat) 1698 { 1699 PetscFunctionBegin; 1700 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1701 if (nump) *nump = ts->num_parameters; 1702 if (Smat) *Smat = ts->mat_sensip; 1703 PetscFunctionReturn(0); 1704 } 1705 1706 /*@ 1707 TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1708 1709 Collective on TS 1710 1711 Input Parameter: 1712 . ts - time stepping context 1713 1714 Level: advanced 1715 1716 Notes: 1717 This function cannot be called until TSStep() has been completed. 1718 1719 .seealso: `TSSolve()`, `TSAdjointCostIntegral()` 1720 @*/ 1721 PetscErrorCode TSForwardCostIntegral(TS ts) 1722 { 1723 PetscFunctionBegin; 1724 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1725 PetscUseTypeMethod(ts, forwardintegral); 1726 PetscFunctionReturn(0); 1727 } 1728 1729 /*@ 1730 TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1731 1732 Collective on TS 1733 1734 Input Parameters: 1735 + ts - the TS context obtained from TSCreate() 1736 - didp - parametric sensitivities of the initial condition 1737 1738 Level: intermediate 1739 1740 Notes: TSSolve() allows users to pass the initial solution directly to TS. But the tangent linear variables cannot be initialized in this way. This function is used to set initial values for tangent linear variables. 1741 1742 .seealso: `TSForwardSetSensitivities()` 1743 @*/ 1744 PetscErrorCode TSForwardSetInitialSensitivities(TS ts, Mat didp) 1745 { 1746 PetscFunctionBegin; 1747 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1748 PetscValidHeaderSpecific(didp, MAT_CLASSID, 2); 1749 if (!ts->mat_sensip) PetscCall(TSForwardSetSensitivities(ts, PETSC_DEFAULT, didp)); 1750 PetscFunctionReturn(0); 1751 } 1752 1753 /*@ 1754 TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 1755 1756 Input Parameter: 1757 . ts - the TS context obtained from TSCreate() 1758 1759 Output Parameters: 1760 + ns - number of stages 1761 - S - tangent linear sensitivities at the intermediate stages 1762 1763 Level: advanced 1764 1765 @*/ 1766 PetscErrorCode TSForwardGetStages(TS ts, PetscInt *ns, Mat **S) 1767 { 1768 PetscFunctionBegin; 1769 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1770 1771 if (!ts->ops->getstages) *S = NULL; 1772 else PetscUseTypeMethod(ts, forwardgetstages, ns, S); 1773 PetscFunctionReturn(0); 1774 } 1775 1776 /*@ 1777 TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 1778 1779 Input Parameters: 1780 + ts - the TS context obtained from TSCreate() 1781 - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1782 1783 Output Parameters: 1784 . quadts - the child TS context 1785 1786 Level: intermediate 1787 1788 .seealso: `TSGetQuadratureTS()` 1789 @*/ 1790 PetscErrorCode TSCreateQuadratureTS(TS ts, PetscBool fwd, TS *quadts) 1791 { 1792 char prefix[128]; 1793 1794 PetscFunctionBegin; 1795 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1796 PetscValidPointer(quadts, 3); 1797 PetscCall(TSDestroy(&ts->quadraturets)); 1798 PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &ts->quadraturets)); 1799 PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets, (PetscObject)ts, 1)); 1800 PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%squad_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "")); 1801 PetscCall(TSSetOptionsPrefix(ts->quadraturets, prefix)); 1802 *quadts = ts->quadraturets; 1803 1804 if (ts->numcost) { 1805 PetscCall(VecCreateSeq(PETSC_COMM_SELF, ts->numcost, &(*quadts)->vec_sol)); 1806 } else { 1807 PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &(*quadts)->vec_sol)); 1808 } 1809 ts->costintegralfwd = fwd; 1810 PetscFunctionReturn(0); 1811 } 1812 1813 /*@ 1814 TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 1815 1816 Input Parameter: 1817 . ts - the TS context obtained from TSCreate() 1818 1819 Output Parameters: 1820 + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1821 - quadts - the child TS context 1822 1823 Level: intermediate 1824 1825 .seealso: `TSCreateQuadratureTS()` 1826 @*/ 1827 PetscErrorCode TSGetQuadratureTS(TS ts, PetscBool *fwd, TS *quadts) 1828 { 1829 PetscFunctionBegin; 1830 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1831 if (fwd) *fwd = ts->costintegralfwd; 1832 if (quadts) *quadts = ts->quadraturets; 1833 PetscFunctionReturn(0); 1834 } 1835 1836 /*@ 1837 TSComputeSNESJacobian - Compute the SNESJacobian 1838 1839 Input Parameters: 1840 + ts - the TS context obtained from TSCreate() 1841 - x - state vector 1842 1843 Output Parameters: 1844 + J - Jacobian matrix 1845 - Jpre - preconditioning matrix for J (may be same as J) 1846 1847 Level: developer 1848 1849 Notes: 1850 Using SNES to compute the Jacobian enables finite differencing when TS Jacobian is not available. 1851 @*/ 1852 PetscErrorCode TSComputeSNESJacobian(TS ts, Vec x, Mat J, Mat Jpre) 1853 { 1854 SNES snes = ts->snes; 1855 PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *) = NULL; 1856 1857 PetscFunctionBegin; 1858 /* 1859 Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object 1860 because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for 1861 explicit methods. Instead, we check the Jacobian compute function directly to determin if FD 1862 coloring is used. 1863 */ 1864 PetscCall(SNESGetJacobian(snes, NULL, NULL, &jac, NULL)); 1865 if (jac == SNESComputeJacobianDefaultColor) { 1866 Vec f; 1867 PetscCall(SNESSetSolution(snes, x)); 1868 PetscCall(SNESGetFunction(snes, &f, NULL, NULL)); 1869 /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */ 1870 PetscCall(SNESComputeFunction(snes, x, f)); 1871 } 1872 PetscCall(SNESComputeJacobian(snes, x, J, Jpre)); 1873 PetscFunctionReturn(0); 1874 } 1875