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