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