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