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