1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/ 2 3 /*@C 4 TaoSetHessian - 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 .seealso: `Tao`, `TaoTypes`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetObjectiveAndGradient()`, `TaoGetHessian()` 28 @*/ 29 PetscErrorCode TaoSetHessian(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 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) tao->user_hessP = ctx; 42 if (func) tao->ops->computehessian = func; 43 if (H) { 44 PetscCall(PetscObjectReference((PetscObject)H)); 45 PetscCall(MatDestroy(&tao->hessian)); 46 tao->hessian = H; 47 } 48 if (Hpre) { 49 PetscCall(PetscObjectReference((PetscObject)Hpre)); 50 PetscCall(MatDestroy(&tao->hessian_pre)); 51 tao->hessian_pre = Hpre; 52 } 53 PetscFunctionReturn(0); 54 } 55 56 /*@C 57 TaoGetHessian - Gets the function to compute the Hessian as well as the location to store the matrix. 58 59 Not collective 60 61 Input Parameter: 62 . tao - the Tao context 63 64 OutputParameters: 65 + H - Matrix used for the hessian 66 . Hpre - Matrix that will be used operated on by preconditioner, can be the same as H 67 . func - Hessian evaluation routine 68 - ctx - user-defined context for private data for the Hessian evaluation routine 69 70 Calling sequence of func: 71 $ func(Tao tao,Vec x,Mat H,Mat Hpre,void *ctx); 72 73 + tao - the Tao context 74 . x - input vector 75 . H - Hessian matrix 76 . Hpre - preconditioner matrix, usually the same as H 77 - ctx - [optional] user-defined Hessian context 78 79 Level: beginner 80 81 .seealso: `Tao`, TaoType`, `TaoGetObjective()`, `TaoGetGradient()`, `TaoGetObjectiveAndGradient()`, `TaoSetHessian()` 82 @*/ 83 PetscErrorCode TaoGetHessian(Tao tao, Mat *H, Mat *Hpre, PetscErrorCode (**func)(Tao, Vec, Mat, Mat, void*), void **ctx) 84 { 85 PetscFunctionBegin; 86 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 87 if (H) *H = tao->hessian; 88 if (Hpre) *Hpre = tao->hessian_pre; 89 if (ctx) *ctx = tao->user_hessP; 90 if (func) *func = tao->ops->computehessian; 91 PetscFunctionReturn(0); 92 } 93 94 PetscErrorCode TaoTestHessian(Tao tao) 95 { 96 Mat A,B,C,D,hessian; 97 Vec x = tao->solution; 98 PetscReal nrm,gnorm; 99 PetscReal threshold = 1.e-5; 100 PetscInt m,n,M,N; 101 PetscBool complete_print = PETSC_FALSE,test = PETSC_FALSE,flg; 102 PetscViewer viewer,mviewer; 103 MPI_Comm comm; 104 PetscInt tabs; 105 static PetscBool directionsprinted = PETSC_FALSE; 106 PetscViewerFormat format; 107 108 PetscFunctionBegin; 109 PetscObjectOptionsBegin((PetscObject)tao); 110 PetscCall(PetscOptionsName("-tao_test_hessian","Compare hand-coded and finite difference Hessians","None",&test)); 111 PetscCall(PetscOptionsReal("-tao_test_hessian", "Threshold for element difference between hand-coded and finite difference being meaningful","None",threshold,&threshold,NULL)); 112 PetscCall(PetscOptionsViewer("-tao_test_hessian_view","View difference between hand-coded and finite difference Hessians element entries","None",&mviewer,&format,&complete_print)); 113 PetscOptionsEnd(); 114 if (!test) PetscFunctionReturn(0); 115 116 PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 117 PetscCall(PetscViewerASCIIGetStdout(comm,&viewer)); 118 PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 119 PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel)); 120 PetscCall(PetscViewerASCIIPrintf(viewer," ---------- Testing Hessian -------------\n")); 121 if (!complete_print && !directionsprinted) { 122 PetscCall(PetscViewerASCIIPrintf(viewer," Run with -tao_test_hessian_view and optionally -tao_test_hessian <threshold> to show difference\n")); 123 PetscCall(PetscViewerASCIIPrintf(viewer," of hand-coded and finite difference Hessian entries greater than <threshold>.\n")); 124 } 125 if (!directionsprinted) { 126 PetscCall(PetscViewerASCIIPrintf(viewer," Testing hand-coded Hessian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n")); 127 PetscCall(PetscViewerASCIIPrintf(viewer," O(1.e-8), the hand-coded Hessian is probably correct.\n")); 128 directionsprinted = PETSC_TRUE; 129 } 130 if (complete_print) PetscCall(PetscViewerPushFormat(mviewer,format)); 131 132 PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian,MATMFFD,&flg)); 133 if (!flg) hessian = tao->hessian; 134 else hessian = tao->hessian_pre; 135 136 while (hessian) { 137 PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)hessian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"")); 138 if (flg) { 139 A = hessian; 140 PetscCall(PetscObjectReference((PetscObject)A)); 141 } else { 142 PetscCall(MatComputeOperator(hessian,MATAIJ,&A)); 143 } 144 145 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); 146 PetscCall(MatGetSize(A,&M,&N)); 147 PetscCall(MatGetLocalSize(A,&m,&n)); 148 PetscCall(MatSetSizes(B,m,n,M,N)); 149 PetscCall(MatSetType(B,((PetscObject)A)->type_name)); 150 PetscCall(MatSetUp(B)); 151 PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 152 153 PetscCall(TaoDefaultComputeHessian(tao,x,B,B,NULL)); 154 155 PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&D)); 156 PetscCall(MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN)); 157 PetscCall(MatNorm(D,NORM_FROBENIUS,&nrm)); 158 PetscCall(MatNorm(A,NORM_FROBENIUS,&gnorm)); 159 PetscCall(MatDestroy(&D)); 160 if (!gnorm) gnorm = 1; /* just in case */ 161 PetscCall(PetscViewerASCIIPrintf(viewer," ||H - Hfd||_F/||H||_F = %g, ||H - Hfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm)); 162 163 if (complete_print) { 164 PetscCall(PetscViewerASCIIPrintf(viewer," Hand-coded Hessian ----------\n")); 165 PetscCall(MatView(A,mviewer)); 166 PetscCall(PetscViewerASCIIPrintf(viewer," Finite difference Hessian ----------\n")); 167 PetscCall(MatView(B,mviewer)); 168 } 169 170 if (complete_print) { 171 PetscInt Istart, Iend, *ccols, bncols, cncols, j, row; 172 PetscScalar *cvals; 173 const PetscInt *bcols; 174 const PetscScalar *bvals; 175 176 PetscCall(MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN)); 177 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 178 PetscCall(MatSetSizes(C,m,n,M,N)); 179 PetscCall(MatSetType(C,((PetscObject)A)->type_name)); 180 PetscCall(MatSetUp(C)); 181 PetscCall(MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 182 PetscCall(MatGetOwnershipRange(B,&Istart,&Iend)); 183 184 for (row = Istart; row < Iend; row++) { 185 PetscCall(MatGetRow(B,row,&bncols,&bcols,&bvals)); 186 PetscCall(PetscMalloc2(bncols,&ccols,bncols,&cvals)); 187 for (j = 0, cncols = 0; j < bncols; j++) { 188 if (PetscAbsScalar(bvals[j]) > threshold) { 189 ccols[cncols] = bcols[j]; 190 cvals[cncols] = bvals[j]; 191 cncols += 1; 192 } 193 } 194 if (cncols) { 195 PetscCall(MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES)); 196 } 197 PetscCall(MatRestoreRow(B,row,&bncols,&bcols,&bvals)); 198 PetscCall(PetscFree2(ccols,cvals)); 199 } 200 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 201 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 202 PetscCall(PetscViewerASCIIPrintf(viewer," Finite-difference minus hand-coded Hessian with tolerance %g ----------\n",(double)threshold)); 203 PetscCall(MatView(C,mviewer)); 204 PetscCall(MatDestroy(&C)); 205 } 206 PetscCall(MatDestroy(&A)); 207 PetscCall(MatDestroy(&B)); 208 209 if (hessian != tao->hessian_pre) { 210 hessian = tao->hessian_pre; 211 PetscCall(PetscViewerASCIIPrintf(viewer," ---------- Testing Hessian for preconditioner -------------\n")); 212 } else hessian = NULL; 213 } 214 if (complete_print) { 215 PetscCall(PetscViewerPopFormat(mviewer)); 216 PetscCall(PetscViewerDestroy(&mviewer)); 217 } 218 PetscCall(PetscViewerASCIISetTab(viewer,tabs)); 219 PetscFunctionReturn(0); 220 } 221 222 /*@C 223 TaoComputeHessian - Computes the Hessian matrix that has been 224 set with `TaoSetHessian()`. 225 226 Collective on tao 227 228 Input Parameters: 229 + tao - the Tao solver context 230 - X - input vector 231 232 Output Parameters: 233 + H - Hessian matrix 234 - Hpre - Preconditioning matrix 235 236 Options Database Keys: 237 + -tao_test_hessian - compare the user provided Hessian with one compute via finite differences to check for errors 238 . -tao_test_hessian <numerical value> - display entries in the difference between the user provided Hessian and finite difference Hessian that are greater than a certain value to help users detect errors 239 - -tao_test_hessian_view - display the user provided Hessian, the finite difference Hessian and the difference between them to help users detect the location of errors in the user provided Hessian 240 241 Notes: 242 Most users should not need to explicitly call this routine, as it 243 is used internally within the minimization solvers. 244 245 `TaoComputeHessian()` is typically used within optimization algorithms, 246 so most users would not generally call this routine 247 themselves. 248 249 Developer Note: 250 The Hessian test mechanism follows `SNESTestJacobian()`. 251 252 Level: developer 253 254 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetHessian()` 255 @*/ 256 PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre) 257 { 258 PetscFunctionBegin; 259 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 260 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 261 PetscCheckSameComm(tao,1,X,2); 262 ++tao->nhess; 263 PetscCall(VecLockReadPush(X)); 264 PetscCall(PetscLogEventBegin(TAO_HessianEval,tao,X,H,Hpre)); 265 PetscCallBack("Tao callback Hessian",(*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP)); 266 PetscCall(PetscLogEventEnd(TAO_HessianEval,tao,X,H,Hpre)); 267 PetscCall(VecLockReadPop(X)); 268 269 PetscCall(TaoTestHessian(tao)); 270 PetscFunctionReturn(0); 271 } 272 273 /*@C 274 TaoComputeJacobian - Computes the Jacobian matrix that has been 275 set with TaoSetJacobianRoutine(). 276 277 Collective on tao 278 279 Input Parameters: 280 + tao - the Tao solver context 281 - X - input vector 282 283 Output Parameters: 284 + J - Jacobian matrix 285 - Jpre - Preconditioning matrix 286 287 Notes: 288 Most users should not need to explicitly call this routine, as it 289 is used internally within the minimization solvers. 290 291 `TaoComputeJacobian()` is typically used within minimization 292 implementations, so most users would not generally call this routine 293 themselves. 294 295 Level: developer 296 297 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianRoutine()` 298 @*/ 299 PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre) 300 { 301 PetscFunctionBegin; 302 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 303 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 304 PetscCheckSameComm(tao,1,X,2); 305 ++tao->njac; 306 PetscCall(VecLockReadPush(X)); 307 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 308 PetscCallBack("Tao callback Jacobian",(*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP)); 309 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 310 PetscCall(VecLockReadPop(X)); 311 PetscFunctionReturn(0); 312 } 313 314 /*@C 315 TaoComputeResidualJacobian - Computes the least-squares residual Jacobian matrix that has been 316 set with `TaoSetJacobianResidual()`. 317 318 Collective on tao 319 320 Input Parameters: 321 + tao - the Tao solver context 322 - X - input vector 323 324 Output Parameters: 325 + J - Jacobian matrix 326 - Jpre - Preconditioning matrix 327 328 Notes: 329 Most users should not need to explicitly call this routine, as it 330 is used internally within the minimization solvers. 331 332 `TaoComputeResidualJacobian()` is typically used within least-squares 333 implementations, so most users would not generally call this routine 334 themselves. 335 336 Level: developer 337 338 .seealso: `Tao`, `TaoComputeResidual()`, `TaoSetJacobianResidual()` 339 @*/ 340 PetscErrorCode TaoComputeResidualJacobian(Tao tao, Vec X, Mat J, Mat Jpre) 341 { 342 PetscFunctionBegin; 343 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 344 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 345 PetscCheckSameComm(tao,1,X,2); 346 ++tao->njac; 347 PetscCall(VecLockReadPush(X)); 348 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 349 PetscCallBack("Tao callback least-squares residual Jacobian",(*tao->ops->computeresidualjacobian)(tao,X,J,Jpre,tao->user_lsjacP)); 350 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 351 PetscCall(VecLockReadPop(X)); 352 PetscFunctionReturn(0); 353 } 354 355 /*@C 356 TaoComputeJacobianState - Computes the Jacobian matrix that has been 357 set with `TaoSetJacobianStateRoutine()`. 358 359 Collective on tao 360 361 Input Parameters: 362 + tao - the Tao solver context 363 - X - input vector 364 365 Output Parameters: 366 + J - Jacobian matrix 367 . Jpre - Preconditioning matrix 368 - Jinv - unknown 369 370 Notes: 371 Most users should not need to explicitly call this routine, as it 372 is used internally within the optimization algorithms. 373 374 Level: developer 375 376 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 377 @*/ 378 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv) 379 { 380 PetscFunctionBegin; 381 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 382 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 383 PetscCheckSameComm(tao,1,X,2); 384 ++tao->njac_state; 385 PetscCall(VecLockReadPush(X)); 386 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 387 PetscCallBack("Tao callback Jacobian(state)",(*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP)); 388 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 389 PetscCall(VecLockReadPop(X)); 390 PetscFunctionReturn(0); 391 } 392 393 /*@C 394 TaoComputeJacobianDesign - Computes the Jacobian matrix that has been 395 set with `TaoSetJacobianDesignRoutine()`. 396 397 Collective on tao 398 399 Input Parameters: 400 + tao - the Tao solver context 401 - X - input vector 402 403 Output Parameters: 404 . J - Jacobian matrix 405 406 Notes: 407 Most users should not need to explicitly call this routine, as it 408 is used internally within the optimization algorithms. 409 410 Level: developer 411 412 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianDesignRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 413 @*/ 414 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J) 415 { 416 PetscFunctionBegin; 417 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 418 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 419 PetscCheckSameComm(tao,1,X,2); 420 ++tao->njac_design; 421 PetscCall(VecLockReadPush(X)); 422 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,NULL)); 423 PetscCallBack("Tao callback Jacobian(design)",(*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP)); 424 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,NULL)); 425 PetscCall(VecLockReadPop(X)); 426 PetscFunctionReturn(0); 427 } 428 429 /*@C 430 TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix. 431 432 Logically collective on tao 433 434 Input Parameters: 435 + tao - the Tao context 436 . J - Matrix used for the jacobian 437 . Jpre - Matrix that will be used operated on by preconditioner, can be same as J 438 . func - Jacobian evaluation routine 439 - ctx - [optional] user-defined context for private data for the 440 Jacobian evaluation routine (may be NULL) 441 442 Calling sequence of func: 443 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 444 445 + tao - the Tao context 446 . x - input vector 447 . J - Jacobian matrix 448 . Jpre - preconditioning matrix, usually the same as J 449 - ctx - [optional] user-defined Jacobian context 450 451 Level: intermediate 452 453 .seealso: `Tao`, `TaoSetGradient()`, `TaoSetObjective()` 454 @*/ 455 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 456 { 457 PetscFunctionBegin; 458 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 459 if (J) { 460 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 461 PetscCheckSameComm(tao,1,J,2); 462 } 463 if (Jpre) { 464 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 465 PetscCheckSameComm(tao,1,Jpre,3); 466 } 467 if (ctx) { 468 tao->user_jacP = ctx; 469 } 470 if (func) { 471 tao->ops->computejacobian = func; 472 } 473 if (J) { 474 PetscCall(PetscObjectReference((PetscObject)J)); 475 PetscCall(MatDestroy(&tao->jacobian)); 476 tao->jacobian = J; 477 } 478 if (Jpre) { 479 PetscCall(PetscObjectReference((PetscObject)Jpre)); 480 PetscCall(MatDestroy(&tao->jacobian_pre)); 481 tao->jacobian_pre=Jpre; 482 } 483 PetscFunctionReturn(0); 484 } 485 486 /*@C 487 TaoSetJacobianResidualRoutine - Sets the function to compute the least-squares residual Jacobian as well as the 488 location to store the matrix. 489 490 Logically collective on tao 491 492 Input Parameters: 493 + tao - the Tao context 494 . J - Matrix used for the jacobian 495 . Jpre - Matrix that will be used operated on by preconditioner, can be same as J 496 . func - Jacobian evaluation routine 497 - ctx - [optional] user-defined context for private data for the 498 Jacobian evaluation routine (may be NULL) 499 500 Calling sequence of func: 501 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 502 503 + tao - the Tao context 504 . x - input vector 505 . J - Jacobian matrix 506 . Jpre - preconditioning matrix, usually the same as J 507 - ctx - [optional] user-defined Jacobian context 508 509 Level: intermediate 510 511 .seealso: `Tao`, `TaoSetGradient()`, `TaoSetObjective()` 512 @*/ 513 PetscErrorCode TaoSetJacobianResidualRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 514 { 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 517 if (J) { 518 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 519 PetscCheckSameComm(tao,1,J,2); 520 } 521 if (Jpre) { 522 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 523 PetscCheckSameComm(tao,1,Jpre,3); 524 } 525 if (ctx) { 526 tao->user_lsjacP = ctx; 527 } 528 if (func) { 529 tao->ops->computeresidualjacobian = func; 530 } 531 if (J) { 532 PetscCall(PetscObjectReference((PetscObject)J)); 533 PetscCall(MatDestroy(&tao->ls_jac)); 534 tao->ls_jac = J; 535 } 536 if (Jpre) { 537 PetscCall(PetscObjectReference((PetscObject)Jpre)); 538 PetscCall(MatDestroy(&tao->ls_jac_pre)); 539 tao->ls_jac_pre=Jpre; 540 } 541 PetscFunctionReturn(0); 542 } 543 544 /*@C 545 TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian 546 (and its inverse) of the constraint function with respect to the state variables. 547 Used only for PDE-constrained optimization. 548 549 Logically collective on tao 550 551 Input Parameters: 552 + tao - the Tao context 553 . J - Matrix used for the jacobian 554 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. Only used if Jinv is NULL 555 . Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse. 556 . func - Jacobian evaluation routine 557 - ctx - [optional] user-defined context for private data for the 558 Jacobian evaluation routine (may be NULL) 559 560 Calling sequence of func: 561 $ func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx); 562 563 + tao - the Tao context 564 . x - input vector 565 . J - Jacobian matrix 566 . Jpre - preconditioner matrix, usually the same as J 567 . Jinv - inverse of J 568 - ctx - [optional] user-defined Jacobian context 569 570 Level: intermediate 571 572 .seealso: `Tao`, `TaoComputeJacobianState()`, `TaoSetJacobianDesignRoutine()`, `TaoSetStateDesignIS()` 573 @*/ 574 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx) 575 { 576 PetscFunctionBegin; 577 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 578 if (J) { 579 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 580 PetscCheckSameComm(tao,1,J,2); 581 } 582 if (Jpre) { 583 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 584 PetscCheckSameComm(tao,1,Jpre,3); 585 } 586 if (Jinv) { 587 PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4); 588 PetscCheckSameComm(tao,1,Jinv,4); 589 } 590 if (ctx) { 591 tao->user_jac_stateP = ctx; 592 } 593 if (func) { 594 tao->ops->computejacobianstate = func; 595 } 596 if (J) { 597 PetscCall(PetscObjectReference((PetscObject)J)); 598 PetscCall(MatDestroy(&tao->jacobian_state)); 599 tao->jacobian_state = J; 600 } 601 if (Jpre) { 602 PetscCall(PetscObjectReference((PetscObject)Jpre)); 603 PetscCall(MatDestroy(&tao->jacobian_state_pre)); 604 tao->jacobian_state_pre=Jpre; 605 } 606 if (Jinv) { 607 PetscCall(PetscObjectReference((PetscObject)Jinv)); 608 PetscCall(MatDestroy(&tao->jacobian_state_inv)); 609 tao->jacobian_state_inv=Jinv; 610 } 611 PetscFunctionReturn(0); 612 } 613 614 /*@C 615 TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of 616 the constraint function with respect to the design variables. Used only for 617 PDE-constrained optimization. 618 619 Logically collective on tao 620 621 Input Parameters: 622 + tao - the Tao context 623 . J - Matrix used for the jacobian 624 . func - Jacobian evaluation routine 625 - ctx - [optional] user-defined context for private data for the 626 Jacobian evaluation routine (may be NULL) 627 628 Calling sequence of func: 629 $ func(Tao tao,Vec x,Mat J,void *ctx); 630 631 + tao - the Tao context 632 . x - input vector 633 . J - Jacobian matrix 634 - ctx - [optional] user-defined Jacobian context 635 636 Level: intermediate 637 638 .seealso: `Tao`, `TaoComputeJacobianDesign()`, `TaoSetJacobianStateRoutine()`, `TaoSetStateDesignIS()` 639 @*/ 640 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx) 641 { 642 PetscFunctionBegin; 643 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 644 if (J) { 645 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 646 PetscCheckSameComm(tao,1,J,2); 647 } 648 if (ctx) { 649 tao->user_jac_designP = ctx; 650 } 651 if (func) { 652 tao->ops->computejacobiandesign = func; 653 } 654 if (J) { 655 PetscCall(PetscObjectReference((PetscObject)J)); 656 PetscCall(MatDestroy(&tao->jacobian_design)); 657 tao->jacobian_design = J; 658 } 659 PetscFunctionReturn(0); 660 } 661 662 /*@ 663 TaoSetStateDesignIS - Indicate to the Tao which variables in the 664 solution vector are state variables and which are design. Only applies to 665 PDE-constrained optimization. 666 667 Logically Collective on Tao 668 669 Input Parameters: 670 + tao - The Tao context 671 . s_is - the index set corresponding to the state variables 672 - d_is - the index set corresponding to the design variables 673 674 Level: intermediate 675 676 .seealso: `Tao`, `TaoSetJacobianStateRoutine()`, `TaoSetJacobianDesignRoutine()` 677 @*/ 678 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is) 679 { 680 PetscFunctionBegin; 681 PetscCall(PetscObjectReference((PetscObject)s_is)); 682 PetscCall(ISDestroy(&tao->state_is)); 683 tao->state_is = s_is; 684 PetscCall(PetscObjectReference((PetscObject)(d_is))); 685 PetscCall(ISDestroy(&tao->design_is)); 686 tao->design_is = d_is; 687 PetscFunctionReturn(0); 688 } 689 690 /*@C 691 TaoComputeJacobianEquality - Computes the Jacobian matrix that has been 692 set with `TaoSetJacobianEqualityRoutine()`. 693 694 Collective on tao 695 696 Input Parameters: 697 + tao - the Tao solver context 698 - X - input vector 699 700 Output Parameters: 701 + J - Jacobian matrix 702 - Jpre - Preconditioning matrix 703 704 Notes: 705 Most users should not need to explicitly call this routine, as it 706 is used internally within the optimization algorithms. 707 708 Level: developer 709 710 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 711 @*/ 712 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre) 713 { 714 PetscFunctionBegin; 715 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 716 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 717 PetscCheckSameComm(tao,1,X,2); 718 ++tao->njac_equality; 719 PetscCall(VecLockReadPush(X)); 720 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 721 PetscCallBack("Tao callback Jacobian(equality)",(*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP)); 722 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 723 PetscCall(VecLockReadPop(X)); 724 PetscFunctionReturn(0); 725 } 726 727 /*@C 728 TaoComputeJacobianInequality - Computes the Jacobian matrix that has been 729 set with `TaoSetJacobianInequalityRoutine()`. 730 731 Collective on tao 732 733 Input Parameters: 734 + tao - the Tao solver context 735 - X - input vector 736 737 Output Parameters: 738 + J - Jacobian matrix 739 - Jpre - Preconditioning matrix 740 741 Notes: 742 Most users should not need to explicitly call this routine, as it 743 is used internally within the minimization solvers. 744 745 Level: developer 746 747 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 748 @*/ 749 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre) 750 { 751 PetscFunctionBegin; 752 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 753 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 754 PetscCheckSameComm(tao,1,X,2); 755 ++tao->njac_inequality; 756 PetscCall(VecLockReadPush(X)); 757 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 758 PetscCallBack("Tao callback Jacobian (inequality)",(*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP)); 759 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 760 PetscCall(VecLockReadPop(X)); 761 PetscFunctionReturn(0); 762 } 763 764 /*@C 765 TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian 766 (and its inverse) of the constraint function with respect to the equality variables. 767 Used only for PDE-constrained optimization. 768 769 Logically collective on tao 770 771 Input Parameters: 772 + tao - the Tao context 773 . J - Matrix used for the jacobian 774 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 775 . func - Jacobian evaluation routine 776 - ctx - [optional] user-defined context for private data for the 777 Jacobian evaluation routine (may be NULL) 778 779 Calling sequence of func: 780 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 781 782 + tao - the Tao context 783 . x - input vector 784 . J - Jacobian matrix 785 . Jpre - preconditioner matrix, usually the same as J 786 - ctx - [optional] user-defined Jacobian context 787 788 Level: intermediate 789 790 .seealso: `Tao`, `TaoComputeJacobianEquality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetEqualityDesignIS()` 791 @*/ 792 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 793 { 794 PetscFunctionBegin; 795 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 796 if (J) { 797 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 798 PetscCheckSameComm(tao,1,J,2); 799 } 800 if (Jpre) { 801 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 802 PetscCheckSameComm(tao,1,Jpre,3); 803 } 804 if (ctx) { 805 tao->user_jac_equalityP = ctx; 806 } 807 if (func) { 808 tao->ops->computejacobianequality = func; 809 } 810 if (J) { 811 PetscCall(PetscObjectReference((PetscObject)J)); 812 PetscCall(MatDestroy(&tao->jacobian_equality)); 813 tao->jacobian_equality = J; 814 } 815 if (Jpre) { 816 PetscCall(PetscObjectReference((PetscObject)Jpre)); 817 PetscCall(MatDestroy(&tao->jacobian_equality_pre)); 818 tao->jacobian_equality_pre=Jpre; 819 } 820 PetscFunctionReturn(0); 821 } 822 823 /*@C 824 TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian 825 (and its inverse) of the constraint function with respect to the inequality variables. 826 Used only for PDE-constrained optimization. 827 828 Logically collective on tao 829 830 Input Parameters: 831 + tao - the Tao context 832 . J - Matrix used for the jacobian 833 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 834 . func - Jacobian evaluation routine 835 - ctx - [optional] user-defined context for private data for the 836 Jacobian evaluation routine (may be NULL) 837 838 Calling sequence of func: 839 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 840 841 + tao - the Tao context 842 . x - input vector 843 . J - Jacobian matrix 844 . Jpre - preconditioner matrix, usually the same as J 845 - ctx - [optional] user-defined Jacobian context 846 847 Level: intermediate 848 849 .seealso: `Tao`, `TaoComputeJacobianInequality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetInequalityDesignIS()` 850 @*/ 851 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx) 852 { 853 PetscFunctionBegin; 854 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 855 if (J) { 856 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 857 PetscCheckSameComm(tao,1,J,2); 858 } 859 if (Jpre) { 860 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 861 PetscCheckSameComm(tao,1,Jpre,3); 862 } 863 if (ctx) { 864 tao->user_jac_inequalityP = ctx; 865 } 866 if (func) { 867 tao->ops->computejacobianinequality = func; 868 } 869 if (J) { 870 PetscCall(PetscObjectReference((PetscObject)J)); 871 PetscCall(MatDestroy(&tao->jacobian_inequality)); 872 tao->jacobian_inequality = J; 873 } 874 if (Jpre) { 875 PetscCall(PetscObjectReference((PetscObject)Jpre)); 876 PetscCall(MatDestroy(&tao->jacobian_inequality_pre)); 877 tao->jacobian_inequality_pre=Jpre; 878 } 879 PetscFunctionReturn(0); 880 } 881