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