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