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