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