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