1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/ 2 3 /*@C 4 TaoSetHessianRoutine - Sets the function to compute the Hessian as well as the location to store the matrix. 5 6 Logically collective on Tao 7 8 Input Parameters: 9 + tao - the Tao context 10 . H - Matrix used for the hessian 11 . Hpre - Matrix that will be used operated on by preconditioner, can be same as H 12 . func - Hessian evaluation routine 13 - ctx - [optional] user-defined context for private data for the 14 Hessian evaluation routine (may be NULL) 15 16 Calling sequence of func: 17 $ func(Tao tao,Vec x,Mat H,Mat Hpre,void *ctx); 18 19 + tao - the Tao context 20 . x - input vector 21 . H - Hessian matrix 22 . Hpre - preconditioner matrix, usually the same as H 23 - ctx - [optional] user-defined Hessian context 24 25 Level: beginner 26 @*/ 27 PetscErrorCode TaoSetHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 28 { 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 33 if (H) { 34 PetscValidHeaderSpecific(H,MAT_CLASSID,2); 35 PetscCheckSameComm(tao,1,H,2); 36 } 37 if (Hpre) { 38 PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 39 PetscCheckSameComm(tao,1,Hpre,3); 40 } 41 if (ctx) { 42 tao->user_hessP = ctx; 43 } 44 if (func) { 45 tao->ops->computehessian = func; 46 } 47 if (H) { 48 ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr); 49 ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr); 50 tao->hessian = H; 51 } 52 if (Hpre) { 53 ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr); 54 ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr); 55 tao->hessian_pre = Hpre; 56 } 57 PetscFunctionReturn(0); 58 } 59 60 /*@C 61 TaoComputeHessian - Computes the Hessian matrix that has been 62 set with TaoSetHessianRoutine(). 63 64 Collective on Tao 65 66 Input Parameters: 67 + tao - the Tao solver context 68 - X - input vector 69 70 Output Parameters: 71 + H - Hessian matrix 72 - Hpre - Preconditioning matrix 73 74 Notes: 75 Most users should not need to explicitly call this routine, as it 76 is used internally within the minimization solvers. 77 78 TaoComputeHessian() is typically used within minimization 79 implementations, so most users would not generally call this routine 80 themselves. 81 82 Level: developer 83 84 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetHessianRoutine() 85 @*/ 86 PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre) 87 { 88 PetscErrorCode ierr; 89 90 PetscFunctionBegin; 91 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 92 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 93 PetscCheckSameComm(tao,1,X,2); 94 if (!tao->ops->computehessian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessianRoutine() first"); 95 ++tao->nhess; 96 ierr = VecLockPush(X);CHKERRQ(ierr); 97 ierr = PetscLogEventBegin(Tao_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr); 98 PetscStackPush("Tao user Hessian function"); 99 ierr = (*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP);CHKERRQ(ierr); 100 PetscStackPop; 101 ierr = PetscLogEventEnd(Tao_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr); 102 ierr = VecLockPop(X);CHKERRQ(ierr); 103 PetscFunctionReturn(0); 104 } 105 106 /*@C 107 TaoComputeJacobian - Computes the Jacobian matrix that has been 108 set with TaoSetJacobianRoutine(). 109 110 Collective on Tao 111 112 Input Parameters: 113 + tao - the Tao solver context 114 - X - input vector 115 116 Output Parameters: 117 + J - Jacobian matrix 118 - Jpre - Preconditioning matrix 119 120 Notes: 121 Most users should not need to explicitly call this routine, as it 122 is used internally within the minimization solvers. 123 124 TaoComputeJacobian() is typically used within minimization 125 implementations, so most users would not generally call this routine 126 themselves. 127 128 Level: developer 129 130 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianRoutine() 131 @*/ 132 PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre) 133 { 134 PetscErrorCode ierr; 135 136 PetscFunctionBegin; 137 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 138 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 139 PetscCheckSameComm(tao,1,X,2); 140 if (!tao->ops->computejacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first"); 141 ++tao->njac; 142 ierr = VecLockPush(X);CHKERRQ(ierr); 143 ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 144 PetscStackPush("Tao user Jacobian function"); 145 ierr = (*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP);CHKERRQ(ierr); 146 PetscStackPop; 147 ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 148 ierr = VecLockPop(X);CHKERRQ(ierr); 149 PetscFunctionReturn(0); 150 } 151 152 /*@C 153 TaoComputeJacobianState - Computes the Jacobian matrix that has been 154 set with TaoSetJacobianStateRoutine(). 155 156 Collective on Tao 157 158 Input Parameters: 159 + tao - the Tao solver context 160 - X - input vector 161 162 Output Parameters: 163 + Jpre - Jacobian matrix 164 - Jinv - Preconditioning matrix 165 166 Notes: 167 Most users should not need to explicitly call this routine, as it 168 is used internally within the minimization solvers. 169 170 TaoComputeJacobianState() is typically used within minimization 171 implementations, so most users would not generally call this routine 172 themselves. 173 174 Level: developer 175 176 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 177 @*/ 178 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv) 179 { 180 PetscErrorCode ierr; 181 182 PetscFunctionBegin; 183 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 184 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 185 PetscCheckSameComm(tao,1,X,2); 186 if (!tao->ops->computejacobianstate) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first"); 187 ++tao->njac_state; 188 ierr = VecLockPush(X);CHKERRQ(ierr); 189 ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 190 PetscStackPush("Tao user Jacobian(state) function"); 191 ierr = (*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP);CHKERRQ(ierr); 192 PetscStackPop; 193 ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 194 ierr = VecLockPop(X);CHKERRQ(ierr); 195 PetscFunctionReturn(0); 196 } 197 198 /*@C 199 TaoComputeJacobianDesign - Computes the Jacobian matrix that has been 200 set with TaoSetJacobianDesignRoutine(). 201 202 Collective on Tao 203 204 Input Parameters: 205 + tao - the Tao solver context 206 - X - input vector 207 208 Output Parameters: 209 . J - Jacobian matrix 210 211 Notes: 212 Most users should not need to explicitly call this routine, as it 213 is used internally within the minimization solvers. 214 215 TaoComputeJacobianDesign() is typically used within minimization 216 implementations, so most users would not generally call this routine 217 themselves. 218 219 Level: developer 220 221 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianDesignRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 222 @*/ 223 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J) 224 { 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 229 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 230 PetscCheckSameComm(tao,1,X,2); 231 if (!tao->ops->computejacobiandesign) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first"); 232 ++tao->njac_design; 233 ierr = VecLockPush(X);CHKERRQ(ierr); 234 ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr); 235 PetscStackPush("Tao user Jacobian(design) function"); 236 ierr = (*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP);CHKERRQ(ierr); 237 PetscStackPop; 238 ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr); 239 ierr = VecLockPop(X);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 /*@C 244 TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix. 245 246 Logically collective on Tao 247 248 Input Parameters: 249 + tao - the Tao context 250 . J - Matrix used for the jacobian 251 . Jpre - Matrix that will be used operated on by preconditioner, can be same as J 252 . func - Jacobian evaluation routine 253 - ctx - [optional] user-defined context for private data for the 254 Jacobian evaluation routine (may be NULL) 255 256 Calling sequence of func: 257 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 258 259 + tao - the Tao context 260 . x - input vector 261 . J - Jacobian matrix 262 . Jpre - preconditioning matrix, usually the same as J 263 - ctx - [optional] user-defined Jacobian context 264 265 Level: intermediate 266 @*/ 267 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 268 { 269 PetscErrorCode ierr; 270 271 PetscFunctionBegin; 272 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 273 if (J) { 274 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 275 PetscCheckSameComm(tao,1,J,2); 276 } 277 if (Jpre) { 278 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 279 PetscCheckSameComm(tao,1,Jpre,3); 280 } 281 if (ctx) { 282 tao->user_jacP = ctx; 283 } 284 if (func) { 285 tao->ops->computejacobian = func; 286 } 287 if (J) { 288 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 289 ierr = MatDestroy(&tao->jacobian);CHKERRQ(ierr); 290 tao->jacobian = J; 291 } 292 if (Jpre) { 293 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 294 ierr = MatDestroy(&tao->jacobian_pre);CHKERRQ(ierr); 295 tao->jacobian_pre=Jpre; 296 } 297 PetscFunctionReturn(0); 298 } 299 300 /*@C 301 TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian 302 (and its inverse) of the constraint function with respect to the state variables. 303 Used only for pde-constrained optimization. 304 305 Logically collective on Tao 306 307 Input Parameters: 308 + tao - the Tao context 309 . J - Matrix used for the jacobian 310 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. Only used if Jinv is NULL 311 . Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse. 312 . func - Jacobian evaluation routine 313 - ctx - [optional] user-defined context for private data for the 314 Jacobian evaluation routine (may be NULL) 315 316 Calling sequence of func: 317 $ func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx); 318 319 + tao - the Tao context 320 . x - input vector 321 . J - Jacobian matrix 322 . Jpre - preconditioner matrix, usually the same as J 323 . Jinv - inverse of J 324 - ctx - [optional] user-defined Jacobian context 325 326 Level: intermediate 327 .seealso: TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS() 328 @*/ 329 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx) 330 { 331 PetscErrorCode ierr; 332 333 PetscFunctionBegin; 334 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 335 if (J) { 336 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 337 PetscCheckSameComm(tao,1,J,2); 338 } 339 if (Jpre) { 340 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 341 PetscCheckSameComm(tao,1,Jpre,3); 342 } 343 if (Jinv) { 344 PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4); 345 PetscCheckSameComm(tao,1,Jinv,4); 346 } 347 if (ctx) { 348 tao->user_jac_stateP = ctx; 349 } 350 if (func) { 351 tao->ops->computejacobianstate = func; 352 } 353 if (J) { 354 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 355 ierr = MatDestroy(&tao->jacobian_state);CHKERRQ(ierr); 356 tao->jacobian_state = J; 357 } 358 if (Jpre) { 359 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 360 ierr = MatDestroy(&tao->jacobian_state_pre);CHKERRQ(ierr); 361 tao->jacobian_state_pre=Jpre; 362 } 363 if (Jinv) { 364 ierr = PetscObjectReference((PetscObject)Jinv);CHKERRQ(ierr); 365 ierr = MatDestroy(&tao->jacobian_state_inv);CHKERRQ(ierr); 366 tao->jacobian_state_inv=Jinv; 367 } 368 PetscFunctionReturn(0); 369 } 370 371 /*@C 372 TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of 373 the constraint function with respect to the design variables. Used only for 374 pde-constrained optimization. 375 376 Logically collective on Tao 377 378 Input Parameters: 379 + tao - the Tao context 380 . J - Matrix used for the jacobian 381 . func - Jacobian evaluation routine 382 - ctx - [optional] user-defined context for private data for the 383 Jacobian evaluation routine (may be NULL) 384 385 Calling sequence of func: 386 $ func(Tao tao,Vec x,Mat J,void *ctx); 387 388 + tao - the Tao context 389 . x - input vector 390 . J - Jacobian matrix 391 - ctx - [optional] user-defined Jacobian context 392 393 Level: intermediate 394 395 .seealso: TaoComputeJacobianDesign(), TaoSetJacobianStateRoutine(), TaoSetStateDesignIS() 396 @*/ 397 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx) 398 { 399 PetscErrorCode ierr; 400 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 403 if (J) { 404 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 405 PetscCheckSameComm(tao,1,J,2); 406 } 407 if (ctx) { 408 tao->user_jac_designP = ctx; 409 } 410 if (func) { 411 tao->ops->computejacobiandesign = func; 412 } 413 if (J) { 414 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 415 ierr = MatDestroy(&tao->jacobian_design);CHKERRQ(ierr); 416 tao->jacobian_design = J; 417 } 418 PetscFunctionReturn(0); 419 } 420 421 /*@ 422 TaoSetStateDesignIS - Indicate to the Tao which variables in the 423 solution vector are state variables and which are design. Only applies to 424 pde-constrained optimization. 425 426 Logically Collective on Tao 427 428 Input Parameters: 429 + tao - The Tao context 430 . s_is - the index set corresponding to the state variables 431 - d_is - the index set corresponding to the design variables 432 433 Level: intermediate 434 435 .seealso: TaoSetJacobianStateRoutine(), TaoSetJacobianDesignRoutine() 436 @*/ 437 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is) 438 { 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 ierr = PetscObjectReference((PetscObject)s_is);CHKERRQ(ierr); 443 ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr); 444 tao->state_is = s_is; 445 ierr = PetscObjectReference((PetscObject)(d_is));CHKERRQ(ierr); 446 ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr); 447 tao->design_is = d_is; 448 PetscFunctionReturn(0); 449 } 450 451 /*@C 452 TaoComputeJacobianEquality - Computes the Jacobian matrix that has been 453 set with TaoSetJacobianEqualityRoutine(). 454 455 Collective on Tao 456 457 Input Parameters: 458 + tao - the Tao solver context 459 - X - input vector 460 461 Output Parameters: 462 + J - Jacobian matrix 463 - Jpre - Preconditioning matrix 464 465 Notes: 466 Most users should not need to explicitly call this routine, as it 467 is used internally within the minimization solvers. 468 469 Level: developer 470 471 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 472 @*/ 473 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre) 474 { 475 PetscErrorCode ierr; 476 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 479 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 480 PetscCheckSameComm(tao,1,X,2); 481 if (!tao->ops->computejacobianequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first"); 482 ++tao->njac_equality; 483 ierr = VecLockPush(X);CHKERRQ(ierr); 484 ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 485 PetscStackPush("Tao user Jacobian(equality) function"); 486 ierr = (*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP);CHKERRQ(ierr); 487 PetscStackPop; 488 ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 489 ierr = VecLockPop(X);CHKERRQ(ierr); 490 PetscFunctionReturn(0); 491 } 492 493 /*@C 494 TaoComputeJacobianInequality - Computes the Jacobian matrix that has been 495 set with TaoSetJacobianInequalityRoutine(). 496 497 Collective on Tao 498 499 Input Parameters: 500 + tao - the Tao solver context 501 - X - input vector 502 503 Output Parameters: 504 + J - Jacobian matrix 505 - Jpre - Preconditioning matrix 506 507 Notes: 508 Most users should not need to explicitly call this routine, as it 509 is used internally within the minimization solvers. 510 511 Level: developer 512 513 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 514 @*/ 515 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre) 516 { 517 PetscErrorCode ierr; 518 519 PetscFunctionBegin; 520 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 521 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 522 PetscCheckSameComm(tao,1,X,2); 523 if (!tao->ops->computejacobianinequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first"); 524 ++tao->njac_inequality; 525 ierr = VecLockPush(X);CHKERRQ(ierr); 526 ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 527 PetscStackPush("Tao user Jacobian(inequality) function"); 528 ierr = (*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP);CHKERRQ(ierr); 529 PetscStackPop; 530 ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 531 ierr = VecLockPop(X);CHKERRQ(ierr); 532 PetscFunctionReturn(0); 533 } 534 535 /*@C 536 TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian 537 (and its inverse) of the constraint function with respect to the equality variables. 538 Used only for pde-constrained optimization. 539 540 Logically collective on Tao 541 542 Input Parameters: 543 + tao - the Tao context 544 . J - Matrix used for the jacobian 545 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 546 . func - Jacobian evaluation routine 547 - ctx - [optional] user-defined context for private data for the 548 Jacobian evaluation routine (may be NULL) 549 550 Calling sequence of func: 551 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 552 553 + tao - the Tao context 554 . x - input vector 555 . J - Jacobian matrix 556 . Jpre - preconditioner matrix, usually the same as J 557 - ctx - [optional] user-defined Jacobian context 558 559 Level: intermediate 560 561 .seealso: TaoComputeJacobianEquality(), TaoSetJacobianDesignRoutine(), TaoSetEqualityDesignIS() 562 @*/ 563 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx) 564 { 565 PetscErrorCode ierr; 566 567 PetscFunctionBegin; 568 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 569 if (J) { 570 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 571 PetscCheckSameComm(tao,1,J,2); 572 } 573 if (Jpre) { 574 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 575 PetscCheckSameComm(tao,1,Jpre,3); 576 } 577 if (ctx) { 578 tao->user_jac_equalityP = ctx; 579 } 580 if (func) { 581 tao->ops->computejacobianequality = func; 582 } 583 if (J) { 584 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 585 ierr = MatDestroy(&tao->jacobian_equality);CHKERRQ(ierr); 586 tao->jacobian_equality = J; 587 } 588 if (Jpre) { 589 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 590 ierr = MatDestroy(&tao->jacobian_equality_pre);CHKERRQ(ierr); 591 tao->jacobian_equality_pre=Jpre; 592 } 593 PetscFunctionReturn(0); 594 } 595 596 /*@C 597 TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian 598 (and its inverse) of the constraint function with respect to the inequality variables. 599 Used only for pde-constrained optimization. 600 601 Logically collective on Tao 602 603 Input Parameters: 604 + tao - the Tao context 605 . J - Matrix used for the jacobian 606 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 607 . func - Jacobian evaluation routine 608 - ctx - [optional] user-defined context for private data for the 609 Jacobian evaluation routine (may be NULL) 610 611 Calling sequence of func: 612 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 613 614 + tao - the Tao context 615 . x - input vector 616 . J - Jacobian matrix 617 . Jpre - preconditioner matrix, usually the same as J 618 - ctx - [optional] user-defined Jacobian context 619 620 Level: intermediate 621 622 .seealso: TaoComputeJacobianInequality(), TaoSetJacobianDesignRoutine(), TaoSetInequalityDesignIS() 623 @*/ 624 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx) 625 { 626 PetscErrorCode ierr; 627 628 PetscFunctionBegin; 629 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 630 if (J) { 631 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 632 PetscCheckSameComm(tao,1,J,2); 633 } 634 if (Jpre) { 635 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 636 PetscCheckSameComm(tao,1,Jpre,3); 637 } 638 if (ctx) { 639 tao->user_jac_inequalityP = ctx; 640 } 641 if (func) { 642 tao->ops->computejacobianinequality = func; 643 } 644 if (J) { 645 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 646 ierr = MatDestroy(&tao->jacobian_inequality);CHKERRQ(ierr); 647 tao->jacobian_inequality = J; 648 } 649 if (Jpre) { 650 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 651 ierr = MatDestroy(&tao->jacobian_inequality_pre);CHKERRQ(ierr); 652 tao->jacobian_inequality_pre=Jpre; 653 } 654 PetscFunctionReturn(0); 655 } 656