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 if (ts->vec_dir) { /* second-order adjoint */ 700 ierr = TSForwardReset(ts);CHKERRQ(ierr); 701 } 702 ts->vecs_sensi = NULL; 703 ts->vecs_sensip = NULL; 704 ts->vecs_sensi2 = NULL; 705 ts->vecs_sensip2 = NULL; 706 ts->vec_dir = NULL; 707 ts->adjointsetupcalled = PETSC_FALSE; 708 PetscFunctionReturn(0); 709 } 710 711 /*@ 712 TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 713 714 Logically Collective on TS 715 716 Input Parameters: 717 + ts - the TS context obtained from TSCreate() 718 . steps - number of steps to use 719 720 Level: intermediate 721 722 Notes: 723 Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 724 so as to integrate back to less than the original timestep 725 726 .keywords: TS, timestep, set, maximum, iterations 727 728 .seealso: TSSetExactFinalTime() 729 @*/ 730 PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 731 { 732 PetscFunctionBegin; 733 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 734 PetscValidLogicalCollectiveInt(ts,steps,2); 735 if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 736 if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 737 ts->adjoint_max_steps = steps; 738 PetscFunctionReturn(0); 739 } 740 741 /*@C 742 TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 743 744 Level: deprecated 745 746 @*/ 747 PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 748 { 749 PetscErrorCode ierr; 750 751 PetscFunctionBegin; 752 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 753 PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 754 755 ts->rhsjacobianp = func; 756 ts->rhsjacobianpctx = ctx; 757 if(Amat) { 758 ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 759 ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 760 ts->Jacp = Amat; 761 } 762 PetscFunctionReturn(0); 763 } 764 765 /*@C 766 TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 767 768 Level: deprecated 769 770 @*/ 771 PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat) 772 { 773 PetscErrorCode ierr; 774 775 PetscFunctionBegin; 776 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 777 PetscValidHeaderSpecific(U,VEC_CLASSID,3); 778 PetscValidPointer(Amat,4); 779 780 PetscStackPush("TS user JacobianP function for sensitivity analysis"); 781 ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 782 PetscStackPop; 783 PetscFunctionReturn(0); 784 } 785 786 /*@ 787 TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction() 788 789 Level: deprecated 790 791 @*/ 792 PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 793 { 794 PetscErrorCode ierr; 795 796 PetscFunctionBegin; 797 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 798 PetscValidHeaderSpecific(U,VEC_CLASSID,3); 799 800 PetscStackPush("TS user DRDY function for sensitivity analysis"); 801 ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 802 PetscStackPop; 803 PetscFunctionReturn(0); 804 } 805 806 /*@ 807 TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction() 808 809 Level: deprecated 810 811 @*/ 812 PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 813 { 814 PetscErrorCode ierr; 815 816 PetscFunctionBegin; 817 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 818 PetscValidHeaderSpecific(U,VEC_CLASSID,3); 819 820 PetscStackPush("TS user DRDP function for sensitivity analysis"); 821 ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 822 PetscStackPop; 823 PetscFunctionReturn(0); 824 } 825 826 /*@C 827 TSAdjointMonitorSensi - monitors the first lambda sensitivity 828 829 Level: intermediate 830 831 .keywords: TS, set, monitor 832 833 .seealso: TSAdjointMonitorSet() 834 @*/ 835 PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 836 { 837 PetscErrorCode ierr; 838 PetscViewer viewer = vf->viewer; 839 840 PetscFunctionBegin; 841 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 842 ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 843 ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 844 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 845 PetscFunctionReturn(0); 846 } 847 848 /*@C 849 TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 850 851 Collective on TS 852 853 Input Parameters: 854 + ts - TS object you wish to monitor 855 . name - the monitor type one is seeking 856 . help - message indicating what monitoring is done 857 . manual - manual page for the monitor 858 . monitor - the monitor function 859 - 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 860 861 Level: developer 862 863 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 864 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 865 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 866 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 867 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 868 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 869 PetscOptionsFList(), PetscOptionsEList() 870 @*/ 871 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*)) 872 { 873 PetscErrorCode ierr; 874 PetscViewer viewer; 875 PetscViewerFormat format; 876 PetscBool flg; 877 878 PetscFunctionBegin; 879 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 880 if (flg) { 881 PetscViewerAndFormat *vf; 882 ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 883 ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 884 if (monitorsetup) { 885 ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 886 } 887 ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 888 } 889 PetscFunctionReturn(0); 890 } 891 892 /*@C 893 TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 894 timestep to display the iteration's progress. 895 896 Logically Collective on TS 897 898 Input Parameters: 899 + ts - the TS context obtained from TSCreate() 900 . adjointmonitor - monitoring routine 901 . adjointmctx - [optional] user-defined context for private data for the 902 monitor routine (use NULL if no context is desired) 903 - adjointmonitordestroy - [optional] routine that frees monitor context 904 (may be NULL) 905 906 Calling sequence of monitor: 907 $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 908 909 + ts - the TS context 910 . 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 911 been interpolated to) 912 . time - current time 913 . u - current iterate 914 . numcost - number of cost functionos 915 . lambda - sensitivities to initial conditions 916 . mu - sensitivities to parameters 917 - adjointmctx - [optional] adjoint monitoring context 918 919 Notes: 920 This routine adds an additional monitor to the list of monitors that 921 already has been loaded. 922 923 Fortran Notes: 924 Only a single monitor function can be set for each TS object 925 926 Level: intermediate 927 928 .keywords: TS, timestep, set, adjoint, monitor 929 930 .seealso: TSAdjointMonitorCancel() 931 @*/ 932 PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 933 { 934 PetscErrorCode ierr; 935 PetscInt i; 936 PetscBool identical; 937 938 PetscFunctionBegin; 939 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 940 for (i=0; i<ts->numbermonitors;i++) { 941 ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 942 if (identical) PetscFunctionReturn(0); 943 } 944 if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 945 ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 946 ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 947 ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 948 PetscFunctionReturn(0); 949 } 950 951 /*@C 952 TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 953 954 Logically Collective on TS 955 956 Input Parameters: 957 . ts - the TS context obtained from TSCreate() 958 959 Notes: 960 There is no way to remove a single, specific monitor. 961 962 Level: intermediate 963 964 .keywords: TS, timestep, set, adjoint, monitor 965 966 .seealso: TSAdjointMonitorSet() 967 @*/ 968 PetscErrorCode TSAdjointMonitorCancel(TS ts) 969 { 970 PetscErrorCode ierr; 971 PetscInt i; 972 973 PetscFunctionBegin; 974 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 975 for (i=0; i<ts->numberadjointmonitors; i++) { 976 if (ts->adjointmonitordestroy[i]) { 977 ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 978 } 979 } 980 ts->numberadjointmonitors = 0; 981 PetscFunctionReturn(0); 982 } 983 984 /*@C 985 TSAdjointMonitorDefault - the default monitor of adjoint computations 986 987 Level: intermediate 988 989 .keywords: TS, set, monitor 990 991 .seealso: TSAdjointMonitorSet() 992 @*/ 993 PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 994 { 995 PetscErrorCode ierr; 996 PetscViewer viewer = vf->viewer; 997 998 PetscFunctionBegin; 999 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1000 ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1001 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1002 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr); 1003 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1004 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1005 PetscFunctionReturn(0); 1006 } 1007 1008 /*@C 1009 TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1010 VecView() for the sensitivities to initial states at each timestep 1011 1012 Collective on TS 1013 1014 Input Parameters: 1015 + ts - the TS context 1016 . step - current time-step 1017 . ptime - current time 1018 . u - current state 1019 . numcost - number of cost functions 1020 . lambda - sensitivities to initial conditions 1021 . mu - sensitivities to parameters 1022 - dummy - either a viewer or NULL 1023 1024 Level: intermediate 1025 1026 .keywords: TS, vector, adjoint, monitor, view 1027 1028 .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 1029 @*/ 1030 PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 1031 { 1032 PetscErrorCode ierr; 1033 TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1034 PetscDraw draw; 1035 PetscReal xl,yl,xr,yr,h; 1036 char time[32]; 1037 1038 PetscFunctionBegin; 1039 if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1040 1041 ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 1042 ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 1043 ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 1044 ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1045 h = yl + .95*(yr - yl); 1046 ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 1047 ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 1048 PetscFunctionReturn(0); 1049 } 1050 1051 /* 1052 TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1053 1054 Collective on TSAdjoint 1055 1056 Input Parameter: 1057 . ts - the TS context 1058 1059 Options Database Keys: 1060 + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1061 . -ts_adjoint_monitor - print information at each adjoint time step 1062 - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1063 1064 Level: developer 1065 1066 Notes: 1067 This is not normally called directly by users 1068 1069 .keywords: TS, trajectory, timestep, set, options, database 1070 1071 .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 1072 */ 1073 PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 1074 { 1075 PetscBool tflg,opt; 1076 PetscErrorCode ierr; 1077 1078 PetscFunctionBegin; 1079 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 1080 ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 1081 tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 1082 ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 1083 if (opt) { 1084 ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 1085 ts->adjoint_solve = tflg; 1086 } 1087 ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 1088 ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 1089 opt = PETSC_FALSE; 1090 ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 1091 if (opt) { 1092 TSMonitorDrawCtx ctx; 1093 PetscInt howoften = 1; 1094 1095 ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 1096 ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 1097 ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 1098 } 1099 PetscFunctionReturn(0); 1100 } 1101 1102 /*@ 1103 TSAdjointStep - Steps one time step backward in the adjoint run 1104 1105 Collective on TS 1106 1107 Input Parameter: 1108 . ts - the TS context obtained from TSCreate() 1109 1110 Level: intermediate 1111 1112 .keywords: TS, adjoint, step 1113 1114 .seealso: TSAdjointSetUp(), TSAdjointSolve() 1115 @*/ 1116 PetscErrorCode TSAdjointStep(TS ts) 1117 { 1118 DM dm; 1119 PetscErrorCode ierr; 1120 1121 PetscFunctionBegin; 1122 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1123 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1124 ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1125 1126 ierr = VecViewFromOptions(ts->vec_sol,(PetscObject)ts,"-ts_view_solution");CHKERRQ(ierr); 1127 1128 ts->reason = TS_CONVERGED_ITERATING; 1129 ts->ptime_prev = ts->ptime; 1130 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); 1131 ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1132 ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 1133 ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1134 ts->adjoint_steps++; ts->steps--; 1135 1136 if (ts->reason < 0) { 1137 if (ts->errorifstepfailed) { 1138 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]); 1139 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]); 1140 else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]); 1141 } 1142 } else if (!ts->reason) { 1143 if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1144 } 1145 PetscFunctionReturn(0); 1146 } 1147 1148 /*@ 1149 TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1150 1151 Collective on TS 1152 1153 Input Parameter: 1154 . ts - the TS context obtained from TSCreate() 1155 1156 Options Database: 1157 . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1158 1159 Level: intermediate 1160 1161 Notes: 1162 This must be called after a call to TSSolve() that solves the forward problem 1163 1164 By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1165 1166 .keywords: TS, timestep, solve 1167 1168 .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 1169 @*/ 1170 PetscErrorCode TSAdjointSolve(TS ts) 1171 { 1172 PetscErrorCode ierr; 1173 1174 PetscFunctionBegin; 1175 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1176 ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1177 1178 /* reset time step and iteration counters */ 1179 ts->adjoint_steps = 0; 1180 ts->ksp_its = 0; 1181 ts->snes_its = 0; 1182 ts->num_snes_failures = 0; 1183 ts->reject = 0; 1184 ts->reason = TS_CONVERGED_ITERATING; 1185 1186 if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1187 if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1188 1189 while (!ts->reason) { 1190 ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1191 ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1192 ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 1193 ierr = TSAdjointStep(ts);CHKERRQ(ierr); 1194 if (ts->vec_costintegral && !ts->costintegralfwd) { 1195 ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 1196 } 1197 } 1198 ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1199 ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1200 ts->solvetime = ts->ptime; 1201 ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 1202 ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 1203 ts->adjoint_max_steps = 0; 1204 PetscFunctionReturn(0); 1205 } 1206 1207 /*@C 1208 TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1209 1210 Collective on TS 1211 1212 Input Parameters: 1213 + ts - time stepping context obtained from TSCreate() 1214 . step - step number that has just completed 1215 . ptime - model time of the state 1216 . u - state at the current model time 1217 . numcost - number of cost functions (dimension of lambda or mu) 1218 . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1219 - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1220 1221 Notes: 1222 TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1223 Users would almost never call this routine directly. 1224 1225 Level: developer 1226 1227 .keywords: TS, timestep 1228 @*/ 1229 PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 1230 { 1231 PetscErrorCode ierr; 1232 PetscInt i,n = ts->numberadjointmonitors; 1233 1234 PetscFunctionBegin; 1235 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1236 PetscValidHeaderSpecific(u,VEC_CLASSID,4); 1237 ierr = VecLockReadPush(u);CHKERRQ(ierr); 1238 for (i=0; i<n; i++) { 1239 ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1240 } 1241 ierr = VecLockReadPop(u);CHKERRQ(ierr); 1242 PetscFunctionReturn(0); 1243 } 1244 1245 /*@ 1246 TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1247 1248 Collective on TS 1249 1250 Input Arguments: 1251 . ts - time stepping context 1252 1253 Level: advanced 1254 1255 Notes: 1256 This function cannot be called until TSAdjointStep() has been completed. 1257 1258 .seealso: TSAdjointSolve(), TSAdjointStep 1259 @*/ 1260 PetscErrorCode TSAdjointCostIntegral(TS ts) 1261 { 1262 PetscErrorCode ierr; 1263 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1264 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); 1265 ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 1266 PetscFunctionReturn(0); 1267 } 1268 1269 /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1270 1271 /*@ 1272 TSForwardSetUp - Sets up the internal data structures for the later use 1273 of forward sensitivity analysis 1274 1275 Collective on TS 1276 1277 Input Parameter: 1278 . ts - the TS context obtained from TSCreate() 1279 1280 Level: advanced 1281 1282 .keywords: TS, forward sensitivity, setup 1283 1284 .seealso: TSCreate(), TSDestroy(), TSSetUp() 1285 @*/ 1286 PetscErrorCode TSForwardSetUp(TS ts) 1287 { 1288 PetscErrorCode ierr; 1289 1290 PetscFunctionBegin; 1291 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1292 if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1293 if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()"); 1294 if (ts->vecs_integral_sensip) { 1295 ierr = VecDuplicateVecs(ts->vec_sol,ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr); 1296 ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr); 1297 } 1298 1299 if (ts->ops->forwardsetup) { 1300 ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 1301 } 1302 ts->forwardsetupcalled = PETSC_TRUE; 1303 PetscFunctionReturn(0); 1304 } 1305 1306 /*@ 1307 TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 1308 1309 Collective on TS 1310 1311 Input Parameter: 1312 . ts - the TS context obtained from TSCreate() 1313 1314 Level: advanced 1315 1316 .keywords: TS, forward sensitivity, reset 1317 1318 .seealso: TSCreate(), TSDestroy(), TSForwardSetUp() 1319 @*/ 1320 PetscErrorCode TSForwardReset(TS ts) 1321 { 1322 PetscErrorCode ierr; 1323 1324 PetscFunctionBegin; 1325 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1326 if (ts->ops->forwardreset) { 1327 ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr); 1328 } 1329 ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1330 ts->vecs_integral_sensip = NULL; 1331 ts->forward_solve = PETSC_FALSE; 1332 ts->forwardsetupcalled = PETSC_FALSE; 1333 PetscFunctionReturn(0); 1334 } 1335 1336 /*@ 1337 TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1338 1339 Input Parameter: 1340 . ts- the TS context obtained from TSCreate() 1341 . numfwdint- number of integrals 1342 . vp = the vectors containing the gradients for each integral w.r.t. parameters 1343 1344 Level: intermediate 1345 1346 .keywords: TS, forward sensitivity 1347 1348 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1349 @*/ 1350 PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 1351 { 1352 PetscFunctionBegin; 1353 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1354 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()"); 1355 if (!ts->numcost) ts->numcost = numfwdint; 1356 1357 ts->vecs_integral_sensip = vp; 1358 PetscFunctionReturn(0); 1359 } 1360 1361 /*@ 1362 TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1363 1364 Input Parameter: 1365 . ts- the TS context obtained from TSCreate() 1366 1367 Output Parameter: 1368 . vp = the vectors containing the gradients for each integral w.r.t. parameters 1369 1370 Level: intermediate 1371 1372 .keywords: TS, forward sensitivity 1373 1374 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1375 @*/ 1376 PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1377 { 1378 PetscFunctionBegin; 1379 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1380 PetscValidPointer(vp,3); 1381 if (numfwdint) *numfwdint = ts->numcost; 1382 if (vp) *vp = ts->vecs_integral_sensip; 1383 PetscFunctionReturn(0); 1384 } 1385 1386 /*@ 1387 TSForwardStep - Compute the forward sensitivity for one time step. 1388 1389 Collective on TS 1390 1391 Input Arguments: 1392 . ts - time stepping context 1393 1394 Level: advanced 1395 1396 Notes: 1397 This function cannot be called until TSStep() has been completed. 1398 1399 .keywords: TS, forward sensitivity 1400 1401 .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1402 @*/ 1403 PetscErrorCode TSForwardStep(TS ts) 1404 { 1405 PetscErrorCode ierr; 1406 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1407 if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1408 ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1409 ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1410 ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1411 PetscFunctionReturn(0); 1412 } 1413 1414 /*@ 1415 TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1416 1417 Logically Collective on TS and Vec 1418 1419 Input Parameters: 1420 + ts - the TS context obtained from TSCreate() 1421 . nump - number of parameters 1422 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1423 1424 Level: beginner 1425 1426 Notes: 1427 Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1428 This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1429 You must call this function before TSSolve(). 1430 The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1431 1432 .keywords: TS, timestep, set, forward sensitivity, initial values 1433 1434 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1435 @*/ 1436 PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1437 { 1438 PetscErrorCode ierr; 1439 1440 PetscFunctionBegin; 1441 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1442 PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1443 ts->forward_solve = PETSC_TRUE; 1444 if (nump == PETSC_DEFAULT) { 1445 ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr); 1446 } else ts->num_parameters = nump; 1447 ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1448 ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1449 ts->mat_sensip = Smat; 1450 PetscFunctionReturn(0); 1451 } 1452 1453 /*@ 1454 TSForwardGetSensitivities - Returns the trajectory sensitivities 1455 1456 Not Collective, but Vec returned is parallel if TS is parallel 1457 1458 Output Parameter: 1459 + ts - the TS context obtained from TSCreate() 1460 . nump - number of parameters 1461 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1462 1463 Level: intermediate 1464 1465 .keywords: TS, forward sensitivity 1466 1467 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1468 @*/ 1469 PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1470 { 1471 PetscFunctionBegin; 1472 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1473 if (nump) *nump = ts->num_parameters; 1474 if (Smat) *Smat = ts->mat_sensip; 1475 PetscFunctionReturn(0); 1476 } 1477 1478 /*@ 1479 TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1480 1481 Collective on TS 1482 1483 Input Arguments: 1484 . ts - time stepping context 1485 1486 Level: advanced 1487 1488 Notes: 1489 This function cannot be called until TSStep() has been completed. 1490 1491 .seealso: TSSolve(), TSAdjointCostIntegral() 1492 @*/ 1493 PetscErrorCode TSForwardCostIntegral(TS ts) 1494 { 1495 PetscErrorCode ierr; 1496 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1497 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); 1498 ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1499 PetscFunctionReturn(0); 1500 } 1501 1502 /*@ 1503 TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1504 1505 Collective on TS and Mat 1506 1507 Input Parameter 1508 + ts - the TS context obtained from TSCreate() 1509 - didp - parametric sensitivities of the initial condition 1510 1511 Level: intermediate 1512 1513 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. 1514 1515 .seealso: TSForwardSetSensitivities() 1516 @*/ 1517 PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp) 1518 { 1519 Vec sp; 1520 PetscInt lsize; 1521 PetscScalar *xarr; 1522 PetscErrorCode ierr; 1523 1524 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1525 if (ts->vec_dir) { /* indicates second-order adjoint caculation */ 1526 Mat A; 1527 ierr = TSForwardGetSensitivities(ts,NULL,&A);CHKERRQ(ierr); 1528 if (!A) { /* create a single-column dense matrix */ 1529 ierr = VecGetLocalSize(ts->vec_dir,&lsize);CHKERRQ(ierr); 1530 ierr = MatCreateDense(PETSC_COMM_WORLD,lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr); 1531 } 1532 ierr = VecDuplicate(ts->vec_dir,&sp);CHKERRQ(ierr); 1533 ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr); 1534 ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr); 1535 if (didp) { 1536 ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr); 1537 } else { /* identity matrix assumed */ 1538 ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr); 1539 } 1540 ierr = VecResetArray(sp);CHKERRQ(ierr); 1541 ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr); 1542 ierr = VecDestroy(&sp);CHKERRQ(ierr); 1543 ierr = TSForwardSetSensitivities(ts,1,A);CHKERRQ(ierr); 1544 } else { 1545 PetscValidHeaderSpecific(didp,MAT_CLASSID,2); 1546 if (!ts->mat_sensip) { 1547 ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr); 1548 } 1549 } 1550 PetscFunctionReturn(0); 1551 } 1552 1553 /*@ 1554 TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 1555 1556 Input Parameter: 1557 . ts - the TS context obtained from TSCreate() 1558 1559 Output Parameters: 1560 + ns - nu 1561 - S - tangent linear sensitivities at the intermediate stages 1562 1563 Level: advanced 1564 1565 .keywords: TS, second-order adjoint, forward sensitivity 1566 @*/ 1567 PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S) 1568 { 1569 PetscErrorCode ierr; 1570 1571 PetscFunctionBegin; 1572 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1573 1574 if (!ts->ops->getstages) *S=NULL; 1575 else { 1576 ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr); 1577 } 1578 PetscFunctionReturn(0); 1579 } 1580