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