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 TSTrajectory tj; 1091 PetscBool match; 1092 PetscErrorCode ierr; 1093 1094 PetscFunctionBegin; 1095 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1096 if (ts->adjointsetupcalled) PetscFunctionReturn(0); 1097 if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 1098 if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 1099 ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr); 1100 ierr = PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match);CHKERRQ(ierr); 1101 if (match) { 1102 PetscBool solution_only; 1103 ierr = TSTrajectoryGetSolutionOnly(tj,&solution_only);CHKERRQ(ierr); 1104 if (solution_only) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"TSAdjoint cannot use the solution-only mode when choosing the Basic TSTrajectory type. Turn it off with -ts_trajectory_solution_only 0"); 1105 } 1106 1107 if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs; 1108 1109 if (ts->quadraturets) { /* if there is integral in the cost function */ 1110 ierr = VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);CHKERRQ(ierr); 1111 if (ts->vecs_sensip){ 1112 ierr = VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr); 1113 } 1114 } 1115 1116 if (ts->ops->adjointsetup) { 1117 ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr); 1118 } 1119 ts->adjointsetupcalled = PETSC_TRUE; 1120 PetscFunctionReturn(0); 1121 } 1122 1123 /*@ 1124 TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 1125 1126 Collective on TS 1127 1128 Input Parameter: 1129 . ts - the TS context obtained from TSCreate() 1130 1131 Level: beginner 1132 1133 .keywords: TS, timestep, reset 1134 1135 .seealso: TSCreate(), TSAdjointSetup(), TSADestroy() 1136 @*/ 1137 PetscErrorCode TSAdjointReset(TS ts) 1138 { 1139 PetscErrorCode ierr; 1140 1141 PetscFunctionBegin; 1142 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1143 if (ts->ops->adjointreset) { 1144 ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr); 1145 } 1146 if (ts->quadraturets) { /* if there is integral in the cost function */ 1147 ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr); 1148 if (ts->vecs_sensip){ 1149 ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr); 1150 } 1151 } 1152 ts->vecs_sensi = NULL; 1153 ts->vecs_sensip = NULL; 1154 ts->vecs_sensi2 = NULL; 1155 ts->vecs_sensi2p = NULL; 1156 ts->vec_dir = NULL; 1157 ts->adjointsetupcalled = PETSC_FALSE; 1158 PetscFunctionReturn(0); 1159 } 1160 1161 /*@ 1162 TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1163 1164 Logically Collective on TS 1165 1166 Input Parameters: 1167 + ts - the TS context obtained from TSCreate() 1168 . steps - number of steps to use 1169 1170 Level: intermediate 1171 1172 Notes: 1173 Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 1174 so as to integrate back to less than the original timestep 1175 1176 .keywords: TS, timestep, set, maximum, iterations 1177 1178 .seealso: TSSetExactFinalTime() 1179 @*/ 1180 PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 1181 { 1182 PetscFunctionBegin; 1183 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1184 PetscValidLogicalCollectiveInt(ts,steps,2); 1185 if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 1186 if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 1187 ts->adjoint_max_steps = steps; 1188 PetscFunctionReturn(0); 1189 } 1190 1191 /*@C 1192 TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1193 1194 Level: deprecated 1195 1196 @*/ 1197 PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 1198 { 1199 PetscErrorCode ierr; 1200 1201 PetscFunctionBegin; 1202 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1203 PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 1204 1205 ts->rhsjacobianp = func; 1206 ts->rhsjacobianpctx = ctx; 1207 if(Amat) { 1208 ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 1209 ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 1210 ts->Jacp = Amat; 1211 } 1212 PetscFunctionReturn(0); 1213 } 1214 1215 /*@C 1216 TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1217 1218 Level: deprecated 1219 1220 @*/ 1221 PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat) 1222 { 1223 PetscErrorCode ierr; 1224 1225 PetscFunctionBegin; 1226 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1227 PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1228 PetscValidPointer(Amat,4); 1229 1230 PetscStackPush("TS user JacobianP function for sensitivity analysis"); 1231 ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 1232 PetscStackPop; 1233 PetscFunctionReturn(0); 1234 } 1235 1236 /*@ 1237 TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction() 1238 1239 Level: deprecated 1240 1241 @*/ 1242 PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 1243 { 1244 PetscErrorCode ierr; 1245 1246 PetscFunctionBegin; 1247 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1248 PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1249 1250 PetscStackPush("TS user DRDY function for sensitivity analysis"); 1251 ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 1252 PetscStackPop; 1253 PetscFunctionReturn(0); 1254 } 1255 1256 /*@ 1257 TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction() 1258 1259 Level: deprecated 1260 1261 @*/ 1262 PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 1263 { 1264 PetscErrorCode ierr; 1265 1266 PetscFunctionBegin; 1267 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1268 PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1269 1270 PetscStackPush("TS user DRDP function for sensitivity analysis"); 1271 ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 1272 PetscStackPop; 1273 PetscFunctionReturn(0); 1274 } 1275 1276 /*@C 1277 TSAdjointMonitorSensi - monitors the first lambda sensitivity 1278 1279 Level: intermediate 1280 1281 .keywords: TS, set, monitor 1282 1283 .seealso: TSAdjointMonitorSet() 1284 @*/ 1285 PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1286 { 1287 PetscErrorCode ierr; 1288 PetscViewer viewer = vf->viewer; 1289 1290 PetscFunctionBegin; 1291 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1292 ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1293 ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 1294 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1295 PetscFunctionReturn(0); 1296 } 1297 1298 /*@C 1299 TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1300 1301 Collective on TS 1302 1303 Input Parameters: 1304 + ts - TS object you wish to monitor 1305 . name - the monitor type one is seeking 1306 . help - message indicating what monitoring is done 1307 . manual - manual page for the monitor 1308 . monitor - the monitor function 1309 - 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 1310 1311 Level: developer 1312 1313 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 1314 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1315 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 1316 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1317 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1318 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1319 PetscOptionsFList(), PetscOptionsEList() 1320 @*/ 1321 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*)) 1322 { 1323 PetscErrorCode ierr; 1324 PetscViewer viewer; 1325 PetscViewerFormat format; 1326 PetscBool flg; 1327 1328 PetscFunctionBegin; 1329 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 1330 if (flg) { 1331 PetscViewerAndFormat *vf; 1332 ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 1333 ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 1334 if (monitorsetup) { 1335 ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 1336 } 1337 ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 1338 } 1339 PetscFunctionReturn(0); 1340 } 1341 1342 /*@C 1343 TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1344 timestep to display the iteration's progress. 1345 1346 Logically Collective on TS 1347 1348 Input Parameters: 1349 + ts - the TS context obtained from TSCreate() 1350 . adjointmonitor - monitoring routine 1351 . adjointmctx - [optional] user-defined context for private data for the 1352 monitor routine (use NULL if no context is desired) 1353 - adjointmonitordestroy - [optional] routine that frees monitor context 1354 (may be NULL) 1355 1356 Calling sequence of monitor: 1357 $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1358 1359 + ts - the TS context 1360 . 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 1361 been interpolated to) 1362 . time - current time 1363 . u - current iterate 1364 . numcost - number of cost functionos 1365 . lambda - sensitivities to initial conditions 1366 . mu - sensitivities to parameters 1367 - adjointmctx - [optional] adjoint monitoring context 1368 1369 Notes: 1370 This routine adds an additional monitor to the list of monitors that 1371 already has been loaded. 1372 1373 Fortran Notes: 1374 Only a single monitor function can be set for each TS object 1375 1376 Level: intermediate 1377 1378 .keywords: TS, timestep, set, adjoint, monitor 1379 1380 .seealso: TSAdjointMonitorCancel() 1381 @*/ 1382 PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 1383 { 1384 PetscErrorCode ierr; 1385 PetscInt i; 1386 PetscBool identical; 1387 1388 PetscFunctionBegin; 1389 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1390 for (i=0; i<ts->numbermonitors;i++) { 1391 ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 1392 if (identical) PetscFunctionReturn(0); 1393 } 1394 if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 1395 ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1396 ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1397 ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 1398 PetscFunctionReturn(0); 1399 } 1400 1401 /*@C 1402 TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1403 1404 Logically Collective on TS 1405 1406 Input Parameters: 1407 . ts - the TS context obtained from TSCreate() 1408 1409 Notes: 1410 There is no way to remove a single, specific monitor. 1411 1412 Level: intermediate 1413 1414 .keywords: TS, timestep, set, adjoint, monitor 1415 1416 .seealso: TSAdjointMonitorSet() 1417 @*/ 1418 PetscErrorCode TSAdjointMonitorCancel(TS ts) 1419 { 1420 PetscErrorCode ierr; 1421 PetscInt i; 1422 1423 PetscFunctionBegin; 1424 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1425 for (i=0; i<ts->numberadjointmonitors; i++) { 1426 if (ts->adjointmonitordestroy[i]) { 1427 ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1428 } 1429 } 1430 ts->numberadjointmonitors = 0; 1431 PetscFunctionReturn(0); 1432 } 1433 1434 /*@C 1435 TSAdjointMonitorDefault - the default monitor of adjoint computations 1436 1437 Level: intermediate 1438 1439 .keywords: TS, set, monitor 1440 1441 .seealso: TSAdjointMonitorSet() 1442 @*/ 1443 PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1444 { 1445 PetscErrorCode ierr; 1446 PetscViewer viewer = vf->viewer; 1447 1448 PetscFunctionBegin; 1449 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1450 ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1451 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1452 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr); 1453 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1454 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1455 PetscFunctionReturn(0); 1456 } 1457 1458 /*@C 1459 TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1460 VecView() for the sensitivities to initial states at each timestep 1461 1462 Collective on TS 1463 1464 Input Parameters: 1465 + ts - the TS context 1466 . step - current time-step 1467 . ptime - current time 1468 . u - current state 1469 . numcost - number of cost functions 1470 . lambda - sensitivities to initial conditions 1471 . mu - sensitivities to parameters 1472 - dummy - either a viewer or NULL 1473 1474 Level: intermediate 1475 1476 .keywords: TS, vector, adjoint, monitor, view 1477 1478 .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 1479 @*/ 1480 PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 1481 { 1482 PetscErrorCode ierr; 1483 TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1484 PetscDraw draw; 1485 PetscReal xl,yl,xr,yr,h; 1486 char time[32]; 1487 1488 PetscFunctionBegin; 1489 if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1490 1491 ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 1492 ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 1493 ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 1494 ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1495 h = yl + .95*(yr - yl); 1496 ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 1497 ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 1498 PetscFunctionReturn(0); 1499 } 1500 1501 /* 1502 TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1503 1504 Collective on TSAdjoint 1505 1506 Input Parameter: 1507 . ts - the TS context 1508 1509 Options Database Keys: 1510 + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1511 . -ts_adjoint_monitor - print information at each adjoint time step 1512 - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1513 1514 Level: developer 1515 1516 Notes: 1517 This is not normally called directly by users 1518 1519 .keywords: TS, trajectory, timestep, set, options, database 1520 1521 .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 1522 */ 1523 PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 1524 { 1525 PetscBool tflg,opt; 1526 PetscErrorCode ierr; 1527 1528 PetscFunctionBegin; 1529 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 1530 ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 1531 tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 1532 ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 1533 if (opt) { 1534 ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 1535 ts->adjoint_solve = tflg; 1536 } 1537 ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 1538 ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 1539 opt = PETSC_FALSE; 1540 ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 1541 if (opt) { 1542 TSMonitorDrawCtx ctx; 1543 PetscInt howoften = 1; 1544 1545 ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 1546 ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 1547 ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 1548 } 1549 PetscFunctionReturn(0); 1550 } 1551 1552 /*@ 1553 TSAdjointStep - Steps one time step backward in the adjoint run 1554 1555 Collective on TS 1556 1557 Input Parameter: 1558 . ts - the TS context obtained from TSCreate() 1559 1560 Level: intermediate 1561 1562 .keywords: TS, adjoint, step 1563 1564 .seealso: TSAdjointSetUp(), TSAdjointSolve() 1565 @*/ 1566 PetscErrorCode TSAdjointStep(TS ts) 1567 { 1568 DM dm; 1569 PetscErrorCode ierr; 1570 1571 PetscFunctionBegin; 1572 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1573 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1574 ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1575 1576 ts->reason = TS_CONVERGED_ITERATING; 1577 ts->ptime_prev = ts->ptime; 1578 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); 1579 ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1580 ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 1581 ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1582 ts->adjoint_steps++; ts->steps--; 1583 1584 if (ts->reason < 0) { 1585 if (ts->errorifstepfailed) { 1586 if (ts->reason == TSADJOINT_DIVERGED_LINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]); 1587 else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]); 1588 } 1589 } else if (!ts->reason) { 1590 if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1591 } 1592 PetscFunctionReturn(0); 1593 } 1594 1595 /*@ 1596 TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1597 1598 Collective on TS 1599 1600 Input Parameter: 1601 . ts - the TS context obtained from TSCreate() 1602 1603 Options Database: 1604 . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1605 1606 Level: intermediate 1607 1608 Notes: 1609 This must be called after a call to TSSolve() that solves the forward problem 1610 1611 By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1612 1613 .keywords: TS, timestep, solve 1614 1615 .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 1616 @*/ 1617 PetscErrorCode TSAdjointSolve(TS ts) 1618 { 1619 PetscErrorCode ierr; 1620 1621 PetscFunctionBegin; 1622 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1623 ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1624 1625 /* reset time step and iteration counters */ 1626 ts->adjoint_steps = 0; 1627 ts->ksp_its = 0; 1628 ts->snes_its = 0; 1629 ts->num_snes_failures = 0; 1630 ts->reject = 0; 1631 ts->reason = TS_CONVERGED_ITERATING; 1632 1633 if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1634 if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1635 1636 while (!ts->reason) { 1637 ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1638 ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1639 ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 1640 ierr = TSAdjointStep(ts);CHKERRQ(ierr); 1641 if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) { 1642 ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 1643 } 1644 } 1645 ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1646 ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1647 ts->solvetime = ts->ptime; 1648 ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 1649 ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 1650 ts->adjoint_max_steps = 0; 1651 PetscFunctionReturn(0); 1652 } 1653 1654 /*@C 1655 TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1656 1657 Collective on TS 1658 1659 Input Parameters: 1660 + ts - time stepping context obtained from TSCreate() 1661 . step - step number that has just completed 1662 . ptime - model time of the state 1663 . u - state at the current model time 1664 . numcost - number of cost functions (dimension of lambda or mu) 1665 . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1666 - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1667 1668 Notes: 1669 TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1670 Users would almost never call this routine directly. 1671 1672 Level: developer 1673 1674 .keywords: TS, timestep 1675 @*/ 1676 PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 1677 { 1678 PetscErrorCode ierr; 1679 PetscInt i,n = ts->numberadjointmonitors; 1680 1681 PetscFunctionBegin; 1682 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1683 PetscValidHeaderSpecific(u,VEC_CLASSID,4); 1684 ierr = VecLockReadPush(u);CHKERRQ(ierr); 1685 for (i=0; i<n; i++) { 1686 ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1687 } 1688 ierr = VecLockReadPop(u);CHKERRQ(ierr); 1689 PetscFunctionReturn(0); 1690 } 1691 1692 /*@ 1693 TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1694 1695 Collective on TS 1696 1697 Input Arguments: 1698 . ts - time stepping context 1699 1700 Level: advanced 1701 1702 Notes: 1703 This function cannot be called until TSAdjointStep() has been completed. 1704 1705 .seealso: TSAdjointSolve(), TSAdjointStep 1706 @*/ 1707 PetscErrorCode TSAdjointCostIntegral(TS ts) 1708 { 1709 PetscErrorCode ierr; 1710 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1711 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); 1712 ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 1713 PetscFunctionReturn(0); 1714 } 1715 1716 /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1717 1718 /*@ 1719 TSForwardSetUp - Sets up the internal data structures for the later use 1720 of forward sensitivity analysis 1721 1722 Collective on TS 1723 1724 Input Parameter: 1725 . ts - the TS context obtained from TSCreate() 1726 1727 Level: advanced 1728 1729 .keywords: TS, forward sensitivity, setup 1730 1731 .seealso: TSCreate(), TSDestroy(), TSSetUp() 1732 @*/ 1733 PetscErrorCode TSForwardSetUp(TS ts) 1734 { 1735 PetscErrorCode ierr; 1736 1737 PetscFunctionBegin; 1738 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1739 if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1740 if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs; 1741 1742 if (ts->ops->forwardsetup) { 1743 ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 1744 } 1745 ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 1746 ts->forwardsetupcalled = PETSC_TRUE; 1747 PetscFunctionReturn(0); 1748 } 1749 1750 /*@ 1751 TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 1752 1753 Collective on TS 1754 1755 Input Parameter: 1756 . ts - the TS context obtained from TSCreate() 1757 1758 Level: advanced 1759 1760 .keywords: TS, forward sensitivity, reset 1761 1762 .seealso: TSCreate(), TSDestroy(), TSForwardSetUp() 1763 @*/ 1764 PetscErrorCode TSForwardReset(TS ts) 1765 { 1766 TS quadts = ts->quadraturets; 1767 PetscErrorCode ierr; 1768 1769 PetscFunctionBegin; 1770 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1771 if (ts->ops->forwardreset) { 1772 ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr); 1773 } 1774 ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1775 if (quadts) { 1776 ierr = MatDestroy(&quadts->mat_sensip);CHKERRQ(ierr); 1777 } 1778 ierr = VecDestroy(&ts->vec_sensip_col);CHKERRQ(ierr); 1779 ts->forward_solve = PETSC_FALSE; 1780 ts->forwardsetupcalled = PETSC_FALSE; 1781 PetscFunctionReturn(0); 1782 } 1783 1784 /*@ 1785 TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1786 1787 Input Parameter: 1788 . ts- the TS context obtained from TSCreate() 1789 . numfwdint- number of integrals 1790 . vp = the vectors containing the gradients for each integral w.r.t. parameters 1791 1792 Level: deprecated 1793 1794 .keywords: TS, forward sensitivity 1795 1796 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1797 @*/ 1798 PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 1799 { 1800 PetscFunctionBegin; 1801 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1802 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()"); 1803 if (!ts->numcost) ts->numcost = numfwdint; 1804 1805 ts->vecs_integral_sensip = vp; 1806 PetscFunctionReturn(0); 1807 } 1808 1809 /*@ 1810 TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1811 1812 Input Parameter: 1813 . ts- the TS context obtained from TSCreate() 1814 1815 Output Parameter: 1816 . vp = the vectors containing the gradients for each integral w.r.t. parameters 1817 1818 Level: deprecated 1819 1820 .keywords: TS, forward sensitivity 1821 1822 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1823 @*/ 1824 PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1825 { 1826 PetscFunctionBegin; 1827 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1828 PetscValidPointer(vp,3); 1829 if (numfwdint) *numfwdint = ts->numcost; 1830 if (vp) *vp = ts->vecs_integral_sensip; 1831 PetscFunctionReturn(0); 1832 } 1833 1834 /*@ 1835 TSForwardStep - Compute the forward sensitivity for one time step. 1836 1837 Collective on TS 1838 1839 Input Arguments: 1840 . ts - time stepping context 1841 1842 Level: advanced 1843 1844 Notes: 1845 This function cannot be called until TSStep() has been completed. 1846 1847 .keywords: TS, forward sensitivity 1848 1849 .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1850 @*/ 1851 PetscErrorCode TSForwardStep(TS ts) 1852 { 1853 PetscErrorCode ierr; 1854 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1855 if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1856 ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1857 ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1858 ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1859 PetscFunctionReturn(0); 1860 } 1861 1862 /*@ 1863 TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1864 1865 Logically Collective on TS and Vec 1866 1867 Input Parameters: 1868 + ts - the TS context obtained from TSCreate() 1869 . nump - number of parameters 1870 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1871 1872 Level: beginner 1873 1874 Notes: 1875 Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1876 This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1877 You must call this function before TSSolve(). 1878 The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1879 1880 .keywords: TS, timestep, set, forward sensitivity, initial values 1881 1882 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1883 @*/ 1884 PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1885 { 1886 PetscErrorCode ierr; 1887 1888 PetscFunctionBegin; 1889 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1890 PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1891 ts->forward_solve = PETSC_TRUE; 1892 if (nump == PETSC_DEFAULT) { 1893 ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr); 1894 } else ts->num_parameters = nump; 1895 ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1896 ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1897 ts->mat_sensip = Smat; 1898 PetscFunctionReturn(0); 1899 } 1900 1901 /*@ 1902 TSForwardGetSensitivities - Returns the trajectory sensitivities 1903 1904 Not Collective, but Vec returned is parallel if TS is parallel 1905 1906 Output Parameter: 1907 + ts - the TS context obtained from TSCreate() 1908 . nump - number of parameters 1909 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1910 1911 Level: intermediate 1912 1913 .keywords: TS, forward sensitivity 1914 1915 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1916 @*/ 1917 PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1918 { 1919 PetscFunctionBegin; 1920 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1921 if (nump) *nump = ts->num_parameters; 1922 if (Smat) *Smat = ts->mat_sensip; 1923 PetscFunctionReturn(0); 1924 } 1925 1926 /*@ 1927 TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1928 1929 Collective on TS 1930 1931 Input Arguments: 1932 . ts - time stepping context 1933 1934 Level: advanced 1935 1936 Notes: 1937 This function cannot be called until TSStep() has been completed. 1938 1939 .seealso: TSSolve(), TSAdjointCostIntegral() 1940 @*/ 1941 PetscErrorCode TSForwardCostIntegral(TS ts) 1942 { 1943 PetscErrorCode ierr; 1944 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1945 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); 1946 ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1947 PetscFunctionReturn(0); 1948 } 1949 1950 /*@ 1951 TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1952 1953 Collective on TS and Mat 1954 1955 Input Parameter 1956 + ts - the TS context obtained from TSCreate() 1957 - didp - parametric sensitivities of the initial condition 1958 1959 Level: intermediate 1960 1961 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. 1962 1963 .seealso: TSForwardSetSensitivities() 1964 @*/ 1965 PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp) 1966 { 1967 PetscErrorCode ierr; 1968 1969 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1970 PetscValidHeaderSpecific(didp,MAT_CLASSID,2); 1971 if (!ts->mat_sensip) { 1972 ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr); 1973 } 1974 PetscFunctionReturn(0); 1975 } 1976 1977 /*@ 1978 TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 1979 1980 Input Parameter: 1981 . ts - the TS context obtained from TSCreate() 1982 1983 Output Parameters: 1984 + ns - number of stages 1985 - S - tangent linear sensitivities at the intermediate stages 1986 1987 Level: advanced 1988 1989 .keywords: TS, second-order adjoint, forward sensitivity 1990 @*/ 1991 PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S) 1992 { 1993 PetscErrorCode ierr; 1994 1995 PetscFunctionBegin; 1996 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1997 1998 if (!ts->ops->getstages) *S=NULL; 1999 else { 2000 ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr); 2001 } 2002 PetscFunctionReturn(0); 2003 } 2004 2005 /*@ 2006 TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 2007 2008 Input Parameter: 2009 + ts - the TS context obtained from TSCreate() 2010 - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 2011 2012 Output Parameters: 2013 . quadts - the child TS context 2014 2015 Level: intermediate 2016 2017 .keywords: TS, quadrature evaluation 2018 2019 .seealso: TSGetQuadratureTS() 2020 @*/ 2021 PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts) 2022 { 2023 char prefix[128]; 2024 PetscErrorCode ierr; 2025 2026 PetscFunctionBegin; 2027 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2028 PetscValidPointer(quadts,2); 2029 ierr = TSDestroy(&ts->quadraturets);CHKERRQ(ierr); 2030 ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);CHKERRQ(ierr); 2031 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);CHKERRQ(ierr); 2032 ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);CHKERRQ(ierr); 2033 ierr = PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");CHKERRQ(ierr); 2034 ierr = TSSetOptionsPrefix(ts->quadraturets,prefix);CHKERRQ(ierr); 2035 *quadts = ts->quadraturets; 2036 2037 if (ts->numcost) { 2038 ierr = VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);CHKERRQ(ierr); 2039 } else { 2040 ierr = VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);CHKERRQ(ierr); 2041 } 2042 ts->costintegralfwd = fwd; 2043 PetscFunctionReturn(0); 2044 } 2045 2046 /*@ 2047 TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 2048 2049 Input Parameter: 2050 . ts - the TS context obtained from TSCreate() 2051 2052 Output Parameters: 2053 + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 2054 - quadts - the child TS context 2055 2056 Level: intermediate 2057 2058 .keywords: TS, quadrature evaluation 2059 2060 .seealso: TSCreateQuadratureTS() 2061 @*/ 2062 PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts) 2063 { 2064 PetscFunctionBegin; 2065 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2066 if (fwd) *fwd = ts->costintegralfwd; 2067 if (quadts) *quadts = ts->quadraturets; 2068 PetscFunctionReturn(0); 2069 } 2070