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