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 PetscCheck(tao->ops->computehessian,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessian() first"); 263 ++tao->nhess; 264 PetscCall(VecLockReadPush(X)); 265 PetscCall(PetscLogEventBegin(TAO_HessianEval,tao,X,H,Hpre)); 266 PetscCallBack("Tao callback Hessian",(*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP)); 267 PetscCall(PetscLogEventEnd(TAO_HessianEval,tao,X,H,Hpre)); 268 PetscCall(VecLockReadPop(X)); 269 270 PetscCall(TaoTestHessian(tao)); 271 PetscFunctionReturn(0); 272 } 273 274 /*@C 275 TaoComputeJacobian - Computes the Jacobian matrix that has been 276 set with TaoSetJacobianRoutine(). 277 278 Collective on tao 279 280 Input Parameters: 281 + tao - the Tao solver context 282 - X - input vector 283 284 Output Parameters: 285 + J - Jacobian matrix 286 - Jpre - Preconditioning matrix 287 288 Notes: 289 Most users should not need to explicitly call this routine, as it 290 is used internally within the minimization solvers. 291 292 `TaoComputeJacobian()` is typically used within minimization 293 implementations, so most users would not generally call this routine 294 themselves. 295 296 Level: developer 297 298 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianRoutine()` 299 @*/ 300 PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre) 301 { 302 PetscFunctionBegin; 303 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 304 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 305 PetscCheckSameComm(tao,1,X,2); 306 PetscCheck(tao->ops->computejacobian,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first"); 307 ++tao->njac; 308 PetscCall(VecLockReadPush(X)); 309 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 310 PetscCallBack("Tao callback Jacobian",(*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP)); 311 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 312 PetscCall(VecLockReadPop(X)); 313 PetscFunctionReturn(0); 314 } 315 316 /*@C 317 TaoComputeResidualJacobian - Computes the least-squares residual Jacobian matrix that has been 318 set with `TaoSetJacobianResidual()`. 319 320 Collective on tao 321 322 Input Parameters: 323 + tao - the Tao solver context 324 - X - input vector 325 326 Output Parameters: 327 + J - Jacobian matrix 328 - Jpre - Preconditioning matrix 329 330 Notes: 331 Most users should not need to explicitly call this routine, as it 332 is used internally within the minimization solvers. 333 334 `TaoComputeResidualJacobian()` is typically used within least-squares 335 implementations, so most users would not generally call this routine 336 themselves. 337 338 Level: developer 339 340 .seealso: `Tao`, `TaoComputeResidual()`, `TaoSetJacobianResidual()` 341 @*/ 342 PetscErrorCode TaoComputeResidualJacobian(Tao tao, Vec X, Mat J, Mat Jpre) 343 { 344 PetscFunctionBegin; 345 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 346 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 347 PetscCheckSameComm(tao,1,X,2); 348 PetscCheck(tao->ops->computeresidualjacobian,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetResidualJacobian() first"); 349 ++tao->njac; 350 PetscCall(VecLockReadPush(X)); 351 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 352 PetscCallBack("Tao callback least-squares residual Jacobian",(*tao->ops->computeresidualjacobian)(tao,X,J,Jpre,tao->user_lsjacP)); 353 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 354 PetscCall(VecLockReadPop(X)); 355 PetscFunctionReturn(0); 356 } 357 358 /*@C 359 TaoComputeJacobianState - Computes the Jacobian matrix that has been 360 set with `TaoSetJacobianStateRoutine()`. 361 362 Collective on tao 363 364 Input Parameters: 365 + tao - the Tao solver context 366 - X - input vector 367 368 Output Parameters: 369 + J - Jacobian matrix 370 . Jpre - Preconditioning matrix 371 - Jinv - unknown 372 373 Notes: 374 Most users should not need to explicitly call this routine, as it 375 is used internally within the optimization algorithms. 376 377 Level: developer 378 379 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 380 @*/ 381 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv) 382 { 383 PetscFunctionBegin; 384 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 385 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 386 PetscCheckSameComm(tao,1,X,2); 387 PetscCheck(tao->ops->computejacobianstate,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first"); 388 ++tao->njac_state; 389 PetscCall(VecLockReadPush(X)); 390 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 391 PetscCallBack("Tao callback Jacobian(state)",(*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP)); 392 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 393 PetscCall(VecLockReadPop(X)); 394 PetscFunctionReturn(0); 395 } 396 397 /*@C 398 TaoComputeJacobianDesign - Computes the Jacobian matrix that has been 399 set with `TaoSetJacobianDesignRoutine()`. 400 401 Collective on tao 402 403 Input Parameters: 404 + tao - the Tao solver context 405 - X - input vector 406 407 Output Parameters: 408 . J - Jacobian matrix 409 410 Notes: 411 Most users should not need to explicitly call this routine, as it 412 is used internally within the optimization algorithms. 413 414 Level: developer 415 416 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianDesignRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 417 @*/ 418 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J) 419 { 420 PetscFunctionBegin; 421 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 422 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 423 PetscCheckSameComm(tao,1,X,2); 424 PetscCheck(tao->ops->computejacobiandesign,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first"); 425 ++tao->njac_design; 426 PetscCall(VecLockReadPush(X)); 427 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,NULL)); 428 PetscCallBack("Tao callback Jacobian(design)",(*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP)); 429 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,NULL)); 430 PetscCall(VecLockReadPop(X)); 431 PetscFunctionReturn(0); 432 } 433 434 /*@C 435 TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix. 436 437 Logically collective on tao 438 439 Input Parameters: 440 + tao - the Tao context 441 . J - Matrix used for the jacobian 442 . Jpre - Matrix that will be used operated on by preconditioner, can be same as J 443 . func - Jacobian evaluation routine 444 - ctx - [optional] user-defined context for private data for the 445 Jacobian evaluation routine (may be NULL) 446 447 Calling sequence of func: 448 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 449 450 + tao - the Tao context 451 . x - input vector 452 . J - Jacobian matrix 453 . Jpre - preconditioning matrix, usually the same as J 454 - ctx - [optional] user-defined Jacobian context 455 456 Level: intermediate 457 458 .seealso: `Tao`, `TaoSetGradient()`, `TaoSetObjective()` 459 @*/ 460 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 461 { 462 PetscFunctionBegin; 463 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 464 if (J) { 465 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 466 PetscCheckSameComm(tao,1,J,2); 467 } 468 if (Jpre) { 469 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 470 PetscCheckSameComm(tao,1,Jpre,3); 471 } 472 if (ctx) { 473 tao->user_jacP = ctx; 474 } 475 if (func) { 476 tao->ops->computejacobian = func; 477 } 478 if (J) { 479 PetscCall(PetscObjectReference((PetscObject)J)); 480 PetscCall(MatDestroy(&tao->jacobian)); 481 tao->jacobian = J; 482 } 483 if (Jpre) { 484 PetscCall(PetscObjectReference((PetscObject)Jpre)); 485 PetscCall(MatDestroy(&tao->jacobian_pre)); 486 tao->jacobian_pre=Jpre; 487 } 488 PetscFunctionReturn(0); 489 } 490 491 /*@C 492 TaoSetJacobianResidualRoutine - Sets the function to compute the least-squares residual Jacobian as well as the 493 location to store the matrix. 494 495 Logically collective on tao 496 497 Input Parameters: 498 + tao - the Tao context 499 . J - Matrix used for the jacobian 500 . Jpre - Matrix that will be used operated on by preconditioner, can be same as J 501 . func - Jacobian evaluation routine 502 - ctx - [optional] user-defined context for private data for the 503 Jacobian evaluation routine (may be NULL) 504 505 Calling sequence of func: 506 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 507 508 + tao - the Tao context 509 . x - input vector 510 . J - Jacobian matrix 511 . Jpre - preconditioning matrix, usually the same as J 512 - ctx - [optional] user-defined Jacobian context 513 514 Level: intermediate 515 516 .seealso: `Tao`, `TaoSetGradient()`, `TaoSetObjective()` 517 @*/ 518 PetscErrorCode TaoSetJacobianResidualRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 519 { 520 PetscFunctionBegin; 521 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 522 if (J) { 523 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 524 PetscCheckSameComm(tao,1,J,2); 525 } 526 if (Jpre) { 527 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 528 PetscCheckSameComm(tao,1,Jpre,3); 529 } 530 if (ctx) { 531 tao->user_lsjacP = ctx; 532 } 533 if (func) { 534 tao->ops->computeresidualjacobian = func; 535 } 536 if (J) { 537 PetscCall(PetscObjectReference((PetscObject)J)); 538 PetscCall(MatDestroy(&tao->ls_jac)); 539 tao->ls_jac = J; 540 } 541 if (Jpre) { 542 PetscCall(PetscObjectReference((PetscObject)Jpre)); 543 PetscCall(MatDestroy(&tao->ls_jac_pre)); 544 tao->ls_jac_pre=Jpre; 545 } 546 PetscFunctionReturn(0); 547 } 548 549 /*@C 550 TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian 551 (and its inverse) of the constraint function with respect to the state variables. 552 Used only for PDE-constrained optimization. 553 554 Logically collective on tao 555 556 Input Parameters: 557 + tao - the Tao context 558 . J - Matrix used for the jacobian 559 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. Only used if Jinv is NULL 560 . Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse. 561 . func - Jacobian evaluation routine 562 - ctx - [optional] user-defined context for private data for the 563 Jacobian evaluation routine (may be NULL) 564 565 Calling sequence of func: 566 $ func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx); 567 568 + tao - the Tao context 569 . x - input vector 570 . J - Jacobian matrix 571 . Jpre - preconditioner matrix, usually the same as J 572 . Jinv - inverse of J 573 - ctx - [optional] user-defined Jacobian context 574 575 Level: intermediate 576 577 .seealso: `Tao`, `TaoComputeJacobianState()`, `TaoSetJacobianDesignRoutine()`, `TaoSetStateDesignIS()` 578 @*/ 579 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx) 580 { 581 PetscFunctionBegin; 582 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 583 if (J) { 584 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 585 PetscCheckSameComm(tao,1,J,2); 586 } 587 if (Jpre) { 588 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 589 PetscCheckSameComm(tao,1,Jpre,3); 590 } 591 if (Jinv) { 592 PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4); 593 PetscCheckSameComm(tao,1,Jinv,4); 594 } 595 if (ctx) { 596 tao->user_jac_stateP = ctx; 597 } 598 if (func) { 599 tao->ops->computejacobianstate = func; 600 } 601 if (J) { 602 PetscCall(PetscObjectReference((PetscObject)J)); 603 PetscCall(MatDestroy(&tao->jacobian_state)); 604 tao->jacobian_state = J; 605 } 606 if (Jpre) { 607 PetscCall(PetscObjectReference((PetscObject)Jpre)); 608 PetscCall(MatDestroy(&tao->jacobian_state_pre)); 609 tao->jacobian_state_pre=Jpre; 610 } 611 if (Jinv) { 612 PetscCall(PetscObjectReference((PetscObject)Jinv)); 613 PetscCall(MatDestroy(&tao->jacobian_state_inv)); 614 tao->jacobian_state_inv=Jinv; 615 } 616 PetscFunctionReturn(0); 617 } 618 619 /*@C 620 TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of 621 the constraint function with respect to the design variables. Used only for 622 PDE-constrained optimization. 623 624 Logically collective on tao 625 626 Input Parameters: 627 + tao - the Tao context 628 . J - Matrix used for the jacobian 629 . func - Jacobian evaluation routine 630 - ctx - [optional] user-defined context for private data for the 631 Jacobian evaluation routine (may be NULL) 632 633 Calling sequence of func: 634 $ func(Tao tao,Vec x,Mat J,void *ctx); 635 636 + tao - the Tao context 637 . x - input vector 638 . J - Jacobian matrix 639 - ctx - [optional] user-defined Jacobian context 640 641 Level: intermediate 642 643 .seealso: `Tao`, `TaoComputeJacobianDesign()`, `TaoSetJacobianStateRoutine()`, `TaoSetStateDesignIS()` 644 @*/ 645 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx) 646 { 647 PetscFunctionBegin; 648 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 649 if (J) { 650 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 651 PetscCheckSameComm(tao,1,J,2); 652 } 653 if (ctx) { 654 tao->user_jac_designP = ctx; 655 } 656 if (func) { 657 tao->ops->computejacobiandesign = func; 658 } 659 if (J) { 660 PetscCall(PetscObjectReference((PetscObject)J)); 661 PetscCall(MatDestroy(&tao->jacobian_design)); 662 tao->jacobian_design = J; 663 } 664 PetscFunctionReturn(0); 665 } 666 667 /*@ 668 TaoSetStateDesignIS - Indicate to the Tao which variables in the 669 solution vector are state variables and which are design. Only applies to 670 PDE-constrained optimization. 671 672 Logically Collective on Tao 673 674 Input Parameters: 675 + tao - The Tao context 676 . s_is - the index set corresponding to the state variables 677 - d_is - the index set corresponding to the design variables 678 679 Level: intermediate 680 681 .seealso: `Tao`, `TaoSetJacobianStateRoutine()`, `TaoSetJacobianDesignRoutine()` 682 @*/ 683 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is) 684 { 685 PetscFunctionBegin; 686 PetscCall(PetscObjectReference((PetscObject)s_is)); 687 PetscCall(ISDestroy(&tao->state_is)); 688 tao->state_is = s_is; 689 PetscCall(PetscObjectReference((PetscObject)(d_is))); 690 PetscCall(ISDestroy(&tao->design_is)); 691 tao->design_is = d_is; 692 PetscFunctionReturn(0); 693 } 694 695 /*@C 696 TaoComputeJacobianEquality - Computes the Jacobian matrix that has been 697 set with `TaoSetJacobianEqualityRoutine()`. 698 699 Collective on tao 700 701 Input Parameters: 702 + tao - the Tao solver context 703 - X - input vector 704 705 Output Parameters: 706 + J - Jacobian matrix 707 - Jpre - Preconditioning matrix 708 709 Notes: 710 Most users should not need to explicitly call this routine, as it 711 is used internally within the optimization algorithms. 712 713 Level: developer 714 715 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 716 @*/ 717 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre) 718 { 719 PetscFunctionBegin; 720 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 721 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 722 PetscCheckSameComm(tao,1,X,2); 723 PetscCheck(tao->ops->computejacobianequality,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first"); 724 ++tao->njac_equality; 725 PetscCall(VecLockReadPush(X)); 726 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 727 PetscCallBack("Tao callback Jacobian(equality)",(*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP)); 728 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 729 PetscCall(VecLockReadPop(X)); 730 PetscFunctionReturn(0); 731 } 732 733 /*@C 734 TaoComputeJacobianInequality - Computes the Jacobian matrix that has been 735 set with `TaoSetJacobianInequalityRoutine()`. 736 737 Collective on tao 738 739 Input Parameters: 740 + tao - the Tao solver context 741 - X - input vector 742 743 Output Parameters: 744 + J - Jacobian matrix 745 - Jpre - Preconditioning matrix 746 747 Notes: 748 Most users should not need to explicitly call this routine, as it 749 is used internally within the minimization solvers. 750 751 Level: developer 752 753 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()` 754 @*/ 755 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre) 756 { 757 PetscFunctionBegin; 758 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 759 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 760 PetscCheckSameComm(tao,1,X,2); 761 PetscCheck(tao->ops->computejacobianinequality,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first"); 762 ++tao->njac_inequality; 763 PetscCall(VecLockReadPush(X)); 764 PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre)); 765 PetscCallBack("Tao callback Jacobian (inequality)",(*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP)); 766 PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre)); 767 PetscCall(VecLockReadPop(X)); 768 PetscFunctionReturn(0); 769 } 770 771 /*@C 772 TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian 773 (and its inverse) of the constraint function with respect to the equality variables. 774 Used only for PDE-constrained optimization. 775 776 Logically collective on tao 777 778 Input Parameters: 779 + tao - the Tao context 780 . J - Matrix used for the jacobian 781 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 782 . func - Jacobian evaluation routine 783 - ctx - [optional] user-defined context for private data for the 784 Jacobian evaluation routine (may be NULL) 785 786 Calling sequence of func: 787 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 788 789 + tao - the Tao context 790 . x - input vector 791 . J - Jacobian matrix 792 . Jpre - preconditioner matrix, usually the same as J 793 - ctx - [optional] user-defined Jacobian context 794 795 Level: intermediate 796 797 .seealso: `Tao`, `TaoComputeJacobianEquality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetEqualityDesignIS()` 798 @*/ 799 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 800 { 801 PetscFunctionBegin; 802 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 803 if (J) { 804 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 805 PetscCheckSameComm(tao,1,J,2); 806 } 807 if (Jpre) { 808 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 809 PetscCheckSameComm(tao,1,Jpre,3); 810 } 811 if (ctx) { 812 tao->user_jac_equalityP = ctx; 813 } 814 if (func) { 815 tao->ops->computejacobianequality = func; 816 } 817 if (J) { 818 PetscCall(PetscObjectReference((PetscObject)J)); 819 PetscCall(MatDestroy(&tao->jacobian_equality)); 820 tao->jacobian_equality = J; 821 } 822 if (Jpre) { 823 PetscCall(PetscObjectReference((PetscObject)Jpre)); 824 PetscCall(MatDestroy(&tao->jacobian_equality_pre)); 825 tao->jacobian_equality_pre=Jpre; 826 } 827 PetscFunctionReturn(0); 828 } 829 830 /*@C 831 TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian 832 (and its inverse) of the constraint function with respect to the inequality variables. 833 Used only for PDE-constrained optimization. 834 835 Logically collective on tao 836 837 Input Parameters: 838 + tao - the Tao context 839 . J - Matrix used for the jacobian 840 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 841 . func - Jacobian evaluation routine 842 - ctx - [optional] user-defined context for private data for the 843 Jacobian evaluation routine (may be NULL) 844 845 Calling sequence of func: 846 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 847 848 + tao - the Tao context 849 . x - input vector 850 . J - Jacobian matrix 851 . Jpre - preconditioner matrix, usually the same as J 852 - ctx - [optional] user-defined Jacobian context 853 854 Level: intermediate 855 856 .seealso: `Tao`, `TaoComputeJacobianInequality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetInequalityDesignIS()` 857 @*/ 858 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx) 859 { 860 PetscFunctionBegin; 861 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 862 if (J) { 863 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 864 PetscCheckSameComm(tao,1,J,2); 865 } 866 if (Jpre) { 867 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 868 PetscCheckSameComm(tao,1,Jpre,3); 869 } 870 if (ctx) { 871 tao->user_jac_inequalityP = ctx; 872 } 873 if (func) { 874 tao->ops->computejacobianinequality = func; 875 } 876 if (J) { 877 PetscCall(PetscObjectReference((PetscObject)J)); 878 PetscCall(MatDestroy(&tao->jacobian_inequality)); 879 tao->jacobian_inequality = J; 880 } 881 if (Jpre) { 882 PetscCall(PetscObjectReference((PetscObject)Jpre)); 883 PetscCall(MatDestroy(&tao->jacobian_inequality_pre)); 884 tao->jacobian_inequality_pre=Jpre; 885 } 886 PetscFunctionReturn(0); 887 } 888