1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/ 2 3 4 /*@C 5 TaoSetHessianRoutine - Sets the function to compute the Hessian as well as the location to store the matrix. 6 7 Logically collective on Tao 8 9 Input Parameters: 10 + tao - the Tao context 11 . H - Matrix used for the hessian 12 . Hpre - Matrix that will be used operated on by preconditioner, can be same as H 13 . func - Hessian evaluation routine 14 - ctx - [optional] user-defined context for private data for the 15 Hessian evaluation routine (may be NULL) 16 17 Calling sequence of func: 18 $ func(Tao tao,Vec x,Mat H,Mat Hpre,void *ctx); 19 20 + tao - the Tao context 21 . x - input vector 22 . H - Hessian matrix 23 . Hpre - preconditioner matrix, usually the same as H 24 - ctx - [optional] user-defined Hessian context 25 26 Level: beginner 27 @*/ 28 PetscErrorCode TaoSetHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 29 { 30 PetscErrorCode ierr; 31 32 PetscFunctionBegin; 33 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 34 if (H) { 35 PetscValidHeaderSpecific(H,MAT_CLASSID,2); 36 PetscCheckSameComm(tao,1,H,2); 37 } 38 if (Hpre) { 39 PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 40 PetscCheckSameComm(tao,1,Hpre,3); 41 } 42 if (ctx) { 43 tao->user_hessP = ctx; 44 } 45 if (func) { 46 tao->ops->computehessian = func; 47 } 48 if (H) { 49 ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr); 50 ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr); 51 tao->hessian = H; 52 } 53 if (Hpre) { 54 ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr); 55 ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr); 56 tao->hessian_pre = Hpre; 57 } 58 PetscFunctionReturn(0); 59 } 60 61 PetscErrorCode TaoTestHessian(Tao tao) 62 { 63 Mat A,B,C,D,hessian; 64 Vec x = tao->solution; 65 PetscErrorCode ierr; 66 PetscReal nrm,gnorm; 67 PetscReal threshold = 1.e-5; 68 PetscInt m,n,M,N; 69 PetscBool complete_print = PETSC_FALSE,test = PETSC_FALSE,flg; 70 PetscViewer viewer,mviewer; 71 MPI_Comm comm; 72 PetscInt tabs; 73 static PetscBool directionsprinted = PETSC_FALSE; 74 PetscViewerFormat format; 75 76 PetscFunctionBegin; 77 ierr = PetscObjectOptionsBegin((PetscObject)tao);CHKERRQ(ierr); 78 ierr = PetscOptionsName("-tao_test_hessian","Compare hand-coded and finite difference Hessians","None",&test);CHKERRQ(ierr); 79 ierr = PetscOptionsReal("-tao_test_hessian", "Threshold for element difference between hand-coded and finite difference being meaningful","None",threshold,&threshold,NULL);CHKERRQ(ierr); 80 ierr = PetscOptionsViewer("-tao_test_hessian_view","View difference between hand-coded and finite difference Hessians element entries","None",&mviewer,&format,&complete_print);CHKERRQ(ierr); 81 ierr = PetscOptionsEnd();CHKERRQ(ierr); 82 if (!test) PetscFunctionReturn(0); 83 84 ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 85 ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 86 ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr); 87 ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel);CHKERRQ(ierr); 88 ierr = PetscViewerASCIIPrintf(viewer," ---------- Testing Hessian -------------\n");CHKERRQ(ierr); 89 if (!complete_print && !directionsprinted) { 90 ierr = PetscViewerASCIIPrintf(viewer," Run with -tao_test_hessian_view and optionally -tao_test_hessian <threshold> to show difference\n");CHKERRQ(ierr); 91 ierr = PetscViewerASCIIPrintf(viewer," of hand-coded and finite difference Hessian entries greater than <threshold>.\n");CHKERRQ(ierr); 92 } 93 if (!directionsprinted) { 94 ierr = PetscViewerASCIIPrintf(viewer," Testing hand-coded Hessian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");CHKERRQ(ierr); 95 ierr = PetscViewerASCIIPrintf(viewer," O(1.e-8), the hand-coded Hessian is probably correct.\n");CHKERRQ(ierr); 96 directionsprinted = PETSC_TRUE; 97 } 98 if (complete_print) { 99 ierr = PetscViewerPushFormat(mviewer,format);CHKERRQ(ierr); 100 } 101 102 ierr = PetscObjectTypeCompare((PetscObject)tao->hessian,MATMFFD,&flg);CHKERRQ(ierr); 103 if (!flg) hessian = tao->hessian; 104 else hessian = tao->hessian_pre; 105 106 while (hessian) { 107 ierr = PetscObjectBaseTypeCompareAny((PetscObject)hessian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr); 108 if (flg) { 109 A = hessian; 110 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 111 } else { 112 ierr = MatComputeOperator(hessian,MATAIJ,&A);CHKERRQ(ierr); 113 } 114 115 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 116 ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 117 ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 118 ierr = MatSetSizes(B,m,n,M,N);CHKERRQ(ierr); 119 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 120 ierr = MatSetUp(B);CHKERRQ(ierr); 121 ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 122 123 ierr = TaoDefaultComputeHessian(tao,x,B,B,NULL);CHKERRQ(ierr); 124 125 ierr = MatDuplicate(B,MAT_COPY_VALUES,&D);CHKERRQ(ierr); 126 ierr = MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 127 ierr = MatNorm(D,NORM_FROBENIUS,&nrm);CHKERRQ(ierr); 128 ierr = MatNorm(A,NORM_FROBENIUS,&gnorm);CHKERRQ(ierr); 129 ierr = MatDestroy(&D);CHKERRQ(ierr); 130 if (!gnorm) gnorm = 1; /* just in case */ 131 ierr = PetscViewerASCIIPrintf(viewer," ||H - Hfd||_F/||H||_F = %g, ||H - Hfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);CHKERRQ(ierr); 132 133 if (complete_print) { 134 ierr = PetscViewerASCIIPrintf(viewer," Hand-coded Hessian ----------\n");CHKERRQ(ierr); 135 ierr = MatView(hessian,mviewer);CHKERRQ(ierr); 136 ierr = PetscViewerASCIIPrintf(viewer," Finite difference Hessian ----------\n");CHKERRQ(ierr); 137 ierr = MatView(B,mviewer);CHKERRQ(ierr); 138 } 139 140 if (complete_print) { 141 PetscInt Istart, Iend, *ccols, bncols, cncols, j, row; 142 PetscScalar *cvals; 143 const PetscInt *bcols; 144 const PetscScalar *bvals; 145 146 ierr = MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 147 ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 148 ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 149 ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 150 ierr = MatSetUp(C);CHKERRQ(ierr); 151 ierr = MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 152 ierr = MatGetOwnershipRange(B,&Istart,&Iend);CHKERRQ(ierr); 153 154 for (row = Istart; row < Iend; row++) { 155 ierr = MatGetRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr); 156 ierr = PetscMalloc2(bncols,&ccols,bncols,&cvals);CHKERRQ(ierr); 157 for (j = 0, cncols = 0; j < bncols; j++) { 158 if (PetscAbsScalar(bvals[j]) > threshold) { 159 ccols[cncols] = bcols[j]; 160 cvals[cncols] = bvals[j]; 161 cncols += 1; 162 } 163 } 164 if (cncols) { 165 ierr = MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);CHKERRQ(ierr); 166 } 167 ierr = MatRestoreRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr); 168 ierr = PetscFree2(ccols,cvals);CHKERRQ(ierr); 169 } 170 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 171 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172 ierr = PetscViewerASCIIPrintf(viewer," Finite-difference minus hand-coded Hessian with tolerance %g ----------\n",(double)threshold);CHKERRQ(ierr); 173 ierr = MatView(C,mviewer);CHKERRQ(ierr); 174 ierr = MatDestroy(&C);CHKERRQ(ierr); 175 } 176 ierr = MatDestroy(&A);CHKERRQ(ierr); 177 ierr = MatDestroy(&B);CHKERRQ(ierr); 178 179 if (hessian != tao->hessian_pre) { 180 hessian = tao->hessian_pre; 181 ierr = PetscViewerASCIIPrintf(viewer," ---------- Testing Hessian for preconditioner -------------\n");CHKERRQ(ierr); 182 } else hessian = NULL; 183 } 184 if (complete_print) { 185 ierr = PetscViewerPopFormat(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 + Jpre - Jacobian matrix 351 - Jinv - Preconditioning matrix 352 353 Notes: 354 Most users should not need to explicitly call this routine, as it 355 is used internally within the minimization solvers. 356 357 TaoComputeJacobianState() is typically used within minimization 358 implementations, so most users would not generally call this routine 359 themselves. 360 361 Level: developer 362 363 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 364 @*/ 365 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv) 366 { 367 PetscErrorCode ierr; 368 369 PetscFunctionBegin; 370 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 371 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 372 PetscCheckSameComm(tao,1,X,2); 373 if (!tao->ops->computejacobianstate) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first"); 374 ++tao->njac_state; 375 ierr = VecLockReadPush(X);CHKERRQ(ierr); 376 ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 377 PetscStackPush("Tao user Jacobian(state) function"); 378 ierr = (*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP);CHKERRQ(ierr); 379 PetscStackPop; 380 ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 381 ierr = VecLockReadPop(X);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 /*@C 386 TaoComputeJacobianDesign - Computes the Jacobian matrix that has been 387 set with TaoSetJacobianDesignRoutine(). 388 389 Collective on Tao 390 391 Input Parameters: 392 + tao - the Tao solver context 393 - X - input vector 394 395 Output Parameters: 396 . J - Jacobian matrix 397 398 Notes: 399 Most users should not need to explicitly call this routine, as it 400 is used internally within the minimization solvers. 401 402 TaoComputeJacobianDesign() is typically used within minimization 403 implementations, so most users would not generally call this routine 404 themselves. 405 406 Level: developer 407 408 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianDesignRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 409 @*/ 410 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J) 411 { 412 PetscErrorCode ierr; 413 414 PetscFunctionBegin; 415 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 416 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 417 PetscCheckSameComm(tao,1,X,2); 418 if (!tao->ops->computejacobiandesign) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first"); 419 ++tao->njac_design; 420 ierr = VecLockReadPush(X);CHKERRQ(ierr); 421 ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr); 422 PetscStackPush("Tao user Jacobian(design) function"); 423 ierr = (*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP);CHKERRQ(ierr); 424 PetscStackPop; 425 ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr); 426 ierr = VecLockReadPop(X);CHKERRQ(ierr); 427 PetscFunctionReturn(0); 428 } 429 430 /*@C 431 TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix. 432 433 Logically collective on Tao 434 435 Input Parameters: 436 + tao - the Tao context 437 . J - Matrix used for the jacobian 438 . Jpre - Matrix that will be used operated on by preconditioner, can be same as J 439 . func - Jacobian evaluation routine 440 - ctx - [optional] user-defined context for private data for the 441 Jacobian evaluation routine (may be NULL) 442 443 Calling sequence of func: 444 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 445 446 + tao - the Tao context 447 . x - input vector 448 . J - Jacobian matrix 449 . Jpre - preconditioning matrix, usually the same as J 450 - ctx - [optional] user-defined Jacobian context 451 452 Level: intermediate 453 @*/ 454 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 455 { 456 PetscErrorCode ierr; 457 458 PetscFunctionBegin; 459 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 460 if (J) { 461 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 462 PetscCheckSameComm(tao,1,J,2); 463 } 464 if (Jpre) { 465 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 466 PetscCheckSameComm(tao,1,Jpre,3); 467 } 468 if (ctx) { 469 tao->user_jacP = ctx; 470 } 471 if (func) { 472 tao->ops->computejacobian = func; 473 } 474 if (J) { 475 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 476 ierr = MatDestroy(&tao->jacobian);CHKERRQ(ierr); 477 tao->jacobian = J; 478 } 479 if (Jpre) { 480 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 481 ierr = MatDestroy(&tao->jacobian_pre);CHKERRQ(ierr); 482 tao->jacobian_pre=Jpre; 483 } 484 PetscFunctionReturn(0); 485 } 486 487 /*@C 488 TaoSetJacobianResidualRoutine - Sets the function to compute the least-squares residual Jacobian as well as the 489 location to store the matrix. 490 491 Logically collective on Tao 492 493 Input Parameters: 494 + tao - the Tao context 495 . J - Matrix used for the jacobian 496 . Jpre - Matrix that will be used operated on by preconditioner, can be same as J 497 . func - Jacobian evaluation routine 498 - ctx - [optional] user-defined context for private data for the 499 Jacobian evaluation routine (may be NULL) 500 501 Calling sequence of func: 502 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 503 504 + tao - the Tao context 505 . x - input vector 506 . J - Jacobian matrix 507 . Jpre - preconditioning matrix, usually the same as J 508 - ctx - [optional] user-defined Jacobian context 509 510 Level: intermediate 511 @*/ 512 PetscErrorCode TaoSetJacobianResidualRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 513 { 514 PetscErrorCode ierr; 515 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 518 if (J) { 519 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 520 PetscCheckSameComm(tao,1,J,2); 521 } 522 if (Jpre) { 523 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 524 PetscCheckSameComm(tao,1,Jpre,3); 525 } 526 if (ctx) { 527 tao->user_lsjacP = ctx; 528 } 529 if (func) { 530 tao->ops->computeresidualjacobian = func; 531 } 532 if (J) { 533 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 534 ierr = MatDestroy(&tao->ls_jac);CHKERRQ(ierr); 535 tao->ls_jac = J; 536 } 537 if (Jpre) { 538 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 539 ierr = MatDestroy(&tao->ls_jac_pre);CHKERRQ(ierr); 540 tao->ls_jac_pre=Jpre; 541 } 542 PetscFunctionReturn(0); 543 } 544 545 /*@C 546 TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian 547 (and its inverse) of the constraint function with respect to the state variables. 548 Used only for pde-constrained optimization. 549 550 Logically collective on Tao 551 552 Input Parameters: 553 + tao - the Tao context 554 . J - Matrix used for the jacobian 555 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. Only used if Jinv is NULL 556 . Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse. 557 . func - Jacobian evaluation routine 558 - ctx - [optional] user-defined context for private data for the 559 Jacobian evaluation routine (may be NULL) 560 561 Calling sequence of func: 562 $ func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx); 563 564 + tao - the Tao context 565 . x - input vector 566 . J - Jacobian matrix 567 . Jpre - preconditioner matrix, usually the same as J 568 . Jinv - inverse of J 569 - ctx - [optional] user-defined Jacobian context 570 571 Level: intermediate 572 .seealso: TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS() 573 @*/ 574 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx) 575 { 576 PetscErrorCode ierr; 577 578 PetscFunctionBegin; 579 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 580 if (J) { 581 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 582 PetscCheckSameComm(tao,1,J,2); 583 } 584 if (Jpre) { 585 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 586 PetscCheckSameComm(tao,1,Jpre,3); 587 } 588 if (Jinv) { 589 PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4); 590 PetscCheckSameComm(tao,1,Jinv,4); 591 } 592 if (ctx) { 593 tao->user_jac_stateP = ctx; 594 } 595 if (func) { 596 tao->ops->computejacobianstate = func; 597 } 598 if (J) { 599 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 600 ierr = MatDestroy(&tao->jacobian_state);CHKERRQ(ierr); 601 tao->jacobian_state = J; 602 } 603 if (Jpre) { 604 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 605 ierr = MatDestroy(&tao->jacobian_state_pre);CHKERRQ(ierr); 606 tao->jacobian_state_pre=Jpre; 607 } 608 if (Jinv) { 609 ierr = PetscObjectReference((PetscObject)Jinv);CHKERRQ(ierr); 610 ierr = MatDestroy(&tao->jacobian_state_inv);CHKERRQ(ierr); 611 tao->jacobian_state_inv=Jinv; 612 } 613 PetscFunctionReturn(0); 614 } 615 616 /*@C 617 TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of 618 the constraint function with respect to the design variables. Used only for 619 pde-constrained optimization. 620 621 Logically collective on Tao 622 623 Input Parameters: 624 + tao - the Tao context 625 . J - Matrix used for the jacobian 626 . func - Jacobian evaluation routine 627 - ctx - [optional] user-defined context for private data for the 628 Jacobian evaluation routine (may be NULL) 629 630 Calling sequence of func: 631 $ func(Tao tao,Vec x,Mat J,void *ctx); 632 633 + tao - the Tao context 634 . x - input vector 635 . J - Jacobian matrix 636 - ctx - [optional] user-defined Jacobian context 637 638 Level: intermediate 639 640 .seealso: TaoComputeJacobianDesign(), TaoSetJacobianStateRoutine(), TaoSetStateDesignIS() 641 @*/ 642 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx) 643 { 644 PetscErrorCode ierr; 645 646 PetscFunctionBegin; 647 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 648 if (J) { 649 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 650 PetscCheckSameComm(tao,1,J,2); 651 } 652 if (ctx) { 653 tao->user_jac_designP = ctx; 654 } 655 if (func) { 656 tao->ops->computejacobiandesign = func; 657 } 658 if (J) { 659 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 660 ierr = MatDestroy(&tao->jacobian_design);CHKERRQ(ierr); 661 tao->jacobian_design = J; 662 } 663 PetscFunctionReturn(0); 664 } 665 666 /*@ 667 TaoSetStateDesignIS - Indicate to the Tao which variables in the 668 solution vector are state variables and which are design. Only applies to 669 pde-constrained optimization. 670 671 Logically Collective on Tao 672 673 Input Parameters: 674 + tao - The Tao context 675 . s_is - the index set corresponding to the state variables 676 - d_is - the index set corresponding to the design variables 677 678 Level: intermediate 679 680 .seealso: TaoSetJacobianStateRoutine(), TaoSetJacobianDesignRoutine() 681 @*/ 682 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is) 683 { 684 PetscErrorCode ierr; 685 686 PetscFunctionBegin; 687 ierr = PetscObjectReference((PetscObject)s_is);CHKERRQ(ierr); 688 ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr); 689 tao->state_is = s_is; 690 ierr = PetscObjectReference((PetscObject)(d_is));CHKERRQ(ierr); 691 ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr); 692 tao->design_is = d_is; 693 PetscFunctionReturn(0); 694 } 695 696 /*@C 697 TaoComputeJacobianEquality - Computes the Jacobian matrix that has been 698 set with TaoSetJacobianEqualityRoutine(). 699 700 Collective on Tao 701 702 Input Parameters: 703 + tao - the Tao solver context 704 - X - input vector 705 706 Output Parameters: 707 + J - Jacobian matrix 708 - Jpre - Preconditioning matrix 709 710 Notes: 711 Most users should not need to explicitly call this routine, as it 712 is used internally within the minimization solvers. 713 714 Level: developer 715 716 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 717 @*/ 718 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre) 719 { 720 PetscErrorCode ierr; 721 722 PetscFunctionBegin; 723 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 724 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 725 PetscCheckSameComm(tao,1,X,2); 726 if (!tao->ops->computejacobianequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first"); 727 ++tao->njac_equality; 728 ierr = VecLockReadPush(X);CHKERRQ(ierr); 729 ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 730 PetscStackPush("Tao user Jacobian(equality) function"); 731 ierr = (*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP);CHKERRQ(ierr); 732 PetscStackPop; 733 ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 734 ierr = VecLockReadPop(X);CHKERRQ(ierr); 735 PetscFunctionReturn(0); 736 } 737 738 /*@C 739 TaoComputeJacobianInequality - Computes the Jacobian matrix that has been 740 set with TaoSetJacobianInequalityRoutine(). 741 742 Collective on Tao 743 744 Input Parameters: 745 + tao - the Tao solver context 746 - X - input vector 747 748 Output Parameters: 749 + J - Jacobian matrix 750 - Jpre - Preconditioning matrix 751 752 Notes: 753 Most users should not need to explicitly call this routine, as it 754 is used internally within the minimization solvers. 755 756 Level: developer 757 758 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS() 759 @*/ 760 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre) 761 { 762 PetscErrorCode ierr; 763 764 PetscFunctionBegin; 765 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 766 PetscValidHeaderSpecific(X, VEC_CLASSID,2); 767 PetscCheckSameComm(tao,1,X,2); 768 if (!tao->ops->computejacobianinequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first"); 769 ++tao->njac_inequality; 770 ierr = VecLockReadPush(X);CHKERRQ(ierr); 771 ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 772 PetscStackPush("Tao user Jacobian(inequality) function"); 773 ierr = (*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP);CHKERRQ(ierr); 774 PetscStackPop; 775 ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr); 776 ierr = VecLockReadPop(X);CHKERRQ(ierr); 777 PetscFunctionReturn(0); 778 } 779 780 /*@C 781 TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian 782 (and its inverse) of the constraint function with respect to the equality variables. 783 Used only for pde-constrained optimization. 784 785 Logically collective on Tao 786 787 Input Parameters: 788 + tao - the Tao context 789 . J - Matrix used for the jacobian 790 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 791 . func - Jacobian evaluation routine 792 - ctx - [optional] user-defined context for private data for the 793 Jacobian evaluation routine (may be NULL) 794 795 Calling sequence of func: 796 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 797 798 + tao - the Tao context 799 . x - input vector 800 . J - Jacobian matrix 801 . Jpre - preconditioner matrix, usually the same as J 802 - ctx - [optional] user-defined Jacobian context 803 804 Level: intermediate 805 806 .seealso: TaoComputeJacobianEquality(), TaoSetJacobianDesignRoutine(), TaoSetEqualityDesignIS() 807 @*/ 808 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx) 809 { 810 PetscErrorCode ierr; 811 812 PetscFunctionBegin; 813 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 814 if (J) { 815 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 816 PetscCheckSameComm(tao,1,J,2); 817 } 818 if (Jpre) { 819 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 820 PetscCheckSameComm(tao,1,Jpre,3); 821 } 822 if (ctx) { 823 tao->user_jac_equalityP = ctx; 824 } 825 if (func) { 826 tao->ops->computejacobianequality = func; 827 } 828 if (J) { 829 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 830 ierr = MatDestroy(&tao->jacobian_equality);CHKERRQ(ierr); 831 tao->jacobian_equality = J; 832 } 833 if (Jpre) { 834 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 835 ierr = MatDestroy(&tao->jacobian_equality_pre);CHKERRQ(ierr); 836 tao->jacobian_equality_pre=Jpre; 837 } 838 PetscFunctionReturn(0); 839 } 840 841 /*@C 842 TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian 843 (and its inverse) of the constraint function with respect to the inequality variables. 844 Used only for pde-constrained optimization. 845 846 Logically collective on Tao 847 848 Input Parameters: 849 + tao - the Tao context 850 . J - Matrix used for the jacobian 851 . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. 852 . func - Jacobian evaluation routine 853 - ctx - [optional] user-defined context for private data for the 854 Jacobian evaluation routine (may be NULL) 855 856 Calling sequence of func: 857 $ func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx); 858 859 + tao - the Tao context 860 . x - input vector 861 . J - Jacobian matrix 862 . Jpre - preconditioner matrix, usually the same as J 863 - ctx - [optional] user-defined Jacobian context 864 865 Level: intermediate 866 867 .seealso: TaoComputeJacobianInequality(), TaoSetJacobianDesignRoutine(), TaoSetInequalityDesignIS() 868 @*/ 869 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx) 870 { 871 PetscErrorCode ierr; 872 873 PetscFunctionBegin; 874 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 875 if (J) { 876 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 877 PetscCheckSameComm(tao,1,J,2); 878 } 879 if (Jpre) { 880 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 881 PetscCheckSameComm(tao,1,Jpre,3); 882 } 883 if (ctx) { 884 tao->user_jac_inequalityP = ctx; 885 } 886 if (func) { 887 tao->ops->computejacobianinequality = func; 888 } 889 if (J) { 890 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 891 ierr = MatDestroy(&tao->jacobian_inequality);CHKERRQ(ierr); 892 tao->jacobian_inequality = J; 893 } 894 if (Jpre) { 895 ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 896 ierr = MatDestroy(&tao->jacobian_inequality_pre);CHKERRQ(ierr); 897 tao->jacobian_inequality_pre=Jpre; 898 } 899 PetscFunctionReturn(0); 900 } 901