1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/ 2 3 /*@ 4 TaoSetSolution - Sets the vector holding the initial guess for the solve 5 6 Logically collective on Tao 7 8 Input Parameters: 9 + tao - the Tao context 10 - x0 - the initial guess 11 12 Level: beginner 13 .seealso: `TaoCreate()`, `TaoSolve()`, `TaoGetSolution()` 14 @*/ 15 PetscErrorCode TaoSetSolution(Tao tao, Vec x0) 16 { 17 PetscFunctionBegin; 18 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 19 if (x0) PetscValidHeaderSpecific(x0,VEC_CLASSID,2); 20 PetscCall(PetscObjectReference((PetscObject)x0)); 21 PetscCall(VecDestroy(&tao->solution)); 22 tao->solution = x0; 23 PetscFunctionReturn(0); 24 } 25 26 PetscErrorCode TaoTestGradient(Tao tao,Vec x,Vec g1) 27 { 28 Vec g2,g3; 29 PetscBool complete_print = PETSC_FALSE,test = PETSC_FALSE; 30 PetscReal hcnorm,fdnorm,hcmax,fdmax,diffmax,diffnorm; 31 PetscScalar dot; 32 MPI_Comm comm; 33 PetscViewer viewer,mviewer; 34 PetscViewerFormat format; 35 PetscInt tabs; 36 static PetscBool directionsprinted = PETSC_FALSE; 37 38 PetscFunctionBegin; 39 PetscObjectOptionsBegin((PetscObject)tao); 40 PetscCall(PetscOptionsName("-tao_test_gradient","Compare hand-coded and finite difference Gradients","None",&test)); 41 PetscCall(PetscOptionsViewer("-tao_test_gradient_view","View difference between hand-coded and finite difference Gradients element entries","None",&mviewer,&format,&complete_print)); 42 PetscOptionsEnd(); 43 if (!test) { 44 if (complete_print) { 45 PetscCall(PetscViewerDestroy(&mviewer)); 46 } 47 PetscFunctionReturn(0); 48 } 49 50 PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 51 PetscCall(PetscViewerASCIIGetStdout(comm,&viewer)); 52 PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 53 PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel)); 54 PetscCall(PetscViewerASCIIPrintf(viewer," ---------- Testing Gradient -------------\n")); 55 if (!complete_print && !directionsprinted) { 56 PetscCall(PetscViewerASCIIPrintf(viewer," Run with -tao_test_gradient_view and optionally -tao_test_gradient <threshold> to show difference\n")); 57 PetscCall(PetscViewerASCIIPrintf(viewer," of hand-coded and finite difference gradient entries greater than <threshold>.\n")); 58 } 59 if (!directionsprinted) { 60 PetscCall(PetscViewerASCIIPrintf(viewer," Testing hand-coded Gradient, if (for double precision runs) ||G - Gfd||/||G|| is\n")); 61 PetscCall(PetscViewerASCIIPrintf(viewer," O(1.e-8), the hand-coded Gradient is probably correct.\n")); 62 directionsprinted = PETSC_TRUE; 63 } 64 if (complete_print) PetscCall(PetscViewerPushFormat(mviewer,format)); 65 66 PetscCall(VecDuplicate(x,&g2)); 67 PetscCall(VecDuplicate(x,&g3)); 68 69 /* Compute finite difference gradient, assume the gradient is already computed by TaoComputeGradient() and put into g1 */ 70 PetscCall(TaoDefaultComputeGradient(tao,x,g2,NULL)); 71 72 PetscCall(VecNorm(g2,NORM_2,&fdnorm)); 73 PetscCall(VecNorm(g1,NORM_2,&hcnorm)); 74 PetscCall(VecNorm(g2,NORM_INFINITY,&fdmax)); 75 PetscCall(VecNorm(g1,NORM_INFINITY,&hcmax)); 76 PetscCall(VecDot(g1,g2,&dot)); 77 PetscCall(VecCopy(g1,g3)); 78 PetscCall(VecAXPY(g3,-1.0,g2)); 79 PetscCall(VecNorm(g3,NORM_2,&diffnorm)); 80 PetscCall(VecNorm(g3,NORM_INFINITY,&diffmax)); 81 PetscCall(PetscViewerASCIIPrintf(viewer," ||Gfd|| %g, ||G|| = %g, angle cosine = (Gfd'G)/||Gfd||||G|| = %g\n", (double)fdnorm, (double)hcnorm, (double)(PetscRealPart(dot)/(fdnorm*hcnorm)))); 82 PetscCall(PetscViewerASCIIPrintf(viewer," 2-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffnorm/PetscMax(hcnorm,fdnorm)),(double)diffnorm)); 83 PetscCall(PetscViewerASCIIPrintf(viewer," max-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffmax/PetscMax(hcmax,fdmax)),(double)diffmax)); 84 85 if (complete_print) { 86 PetscCall(PetscViewerASCIIPrintf(viewer," Hand-coded gradient ----------\n")); 87 PetscCall(VecView(g1,mviewer)); 88 PetscCall(PetscViewerASCIIPrintf(viewer," Finite difference gradient ----------\n")); 89 PetscCall(VecView(g2,mviewer)); 90 PetscCall(PetscViewerASCIIPrintf(viewer," Hand-coded minus finite-difference gradient ----------\n")); 91 PetscCall(VecView(g3,mviewer)); 92 } 93 PetscCall(VecDestroy(&g2)); 94 PetscCall(VecDestroy(&g3)); 95 96 if (complete_print) { 97 PetscCall(PetscViewerPopFormat(mviewer)); 98 PetscCall(PetscViewerDestroy(&mviewer)); 99 } 100 PetscCall(PetscViewerASCIISetTab(viewer,tabs)); 101 PetscFunctionReturn(0); 102 } 103 104 /*@ 105 TaoComputeGradient - Computes the gradient of the objective function 106 107 Collective on Tao 108 109 Input Parameters: 110 + tao - the Tao context 111 - X - input vector 112 113 Output Parameter: 114 . G - gradient vector 115 116 Options Database Keys: 117 + -tao_test_gradient - compare the user provided gradient with one compute via finite differences to check for errors 118 - -tao_test_gradient_view - display the user provided gradient, the finite difference gradient and the difference between them to help users detect the location of errors in the user provided gradient 119 120 Notes: 121 TaoComputeGradient() is typically used within minimization implementations, 122 so most users would not generally call this routine themselves. 123 124 Level: advanced 125 126 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetGradient()` 127 @*/ 128 PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G) 129 { 130 PetscReal dummy; 131 132 PetscFunctionBegin; 133 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 134 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 135 PetscValidHeaderSpecific(G,VEC_CLASSID,3); 136 PetscCheckSameComm(tao,1,X,2); 137 PetscCheckSameComm(tao,1,G,3); 138 PetscCall(VecLockReadPush(X)); 139 if (tao->ops->computegradient) { 140 PetscCall(PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL)); 141 PetscStackPush("Tao user gradient evaluation routine"); 142 PetscCall((*tao->ops->computegradient)(tao,X,G,tao->user_gradP)); 143 PetscStackPop; 144 PetscCall(PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL)); 145 tao->ngrads++; 146 } else if (tao->ops->computeobjectiveandgradient) { 147 PetscCall(PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL)); 148 PetscStackPush("Tao user objective/gradient evaluation routine"); 149 PetscCall((*tao->ops->computeobjectiveandgradient)(tao,X,&dummy,G,tao->user_objgradP)); 150 PetscStackPop; 151 PetscCall(PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL)); 152 tao->nfuncgrads++; 153 } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetGradient() has not been called"); 154 PetscCall(VecLockReadPop(X)); 155 156 PetscCall(TaoTestGradient(tao,X,G)); 157 PetscFunctionReturn(0); 158 } 159 160 /*@ 161 TaoComputeObjective - Computes the objective function value at a given point 162 163 Collective on Tao 164 165 Input Parameters: 166 + tao - the Tao context 167 - X - input vector 168 169 Output Parameter: 170 . f - Objective value at X 171 172 Notes: 173 TaoComputeObjective() is typically used within minimization implementations, 174 so most users would not generally call this routine themselves. 175 176 Level: advanced 177 178 .seealso: `TaoComputeGradient()`, `TaoComputeObjectiveAndGradient()`, `TaoSetObjective()` 179 @*/ 180 PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f) 181 { 182 Vec temp; 183 184 PetscFunctionBegin; 185 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 186 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 187 PetscCheckSameComm(tao,1,X,2); 188 PetscCall(VecLockReadPush(X)); 189 if (tao->ops->computeobjective) { 190 PetscCall(PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL)); 191 PetscStackPush("Tao user objective evaluation routine"); 192 PetscCall((*tao->ops->computeobjective)(tao,X,f,tao->user_objP)); 193 PetscStackPop; 194 PetscCall(PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL)); 195 tao->nfuncs++; 196 } else if (tao->ops->computeobjectiveandgradient) { 197 PetscCall(PetscInfo(tao,"Duplicating variable vector in order to call func/grad routine\n")); 198 PetscCall(VecDuplicate(X,&temp)); 199 PetscCall(PetscLogEventBegin(TAO_ObjGradEval,tao,X,NULL,NULL)); 200 PetscStackPush("Tao user objective/gradient evaluation routine"); 201 PetscCall((*tao->ops->computeobjectiveandgradient)(tao,X,f,temp,tao->user_objgradP)); 202 PetscStackPop; 203 PetscCall(PetscLogEventEnd(TAO_ObjGradEval,tao,X,NULL,NULL)); 204 PetscCall(VecDestroy(&temp)); 205 tao->nfuncgrads++; 206 } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjective() has not been called"); 207 PetscCall(PetscInfo(tao,"TAO Function evaluation: %20.19e\n",(double)(*f))); 208 PetscCall(VecLockReadPop(X)); 209 PetscFunctionReturn(0); 210 } 211 212 /*@ 213 TaoComputeObjectiveAndGradient - Computes the objective function value at a given point 214 215 Collective on Tao 216 217 Input Parameters: 218 + tao - the Tao context 219 - X - input vector 220 221 Output Parameters: 222 + f - Objective value at X 223 - g - Gradient vector at X 224 225 Notes: 226 TaoComputeObjectiveAndGradient() is typically used within minimization implementations, 227 so most users would not generally call this routine themselves. 228 229 Level: advanced 230 231 .seealso: `TaoComputeGradient()`, `TaoComputeObjectiveAndGradient()`, `TaoSetObjective()` 232 @*/ 233 PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G) 234 { 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 237 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 238 PetscValidHeaderSpecific(G,VEC_CLASSID,4); 239 PetscCheckSameComm(tao,1,X,2); 240 PetscCheckSameComm(tao,1,G,4); 241 PetscCall(VecLockReadPush(X)); 242 if (tao->ops->computeobjectiveandgradient) { 243 PetscCall(PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL)); 244 if (tao->ops->computegradient == TaoDefaultComputeGradient) { 245 PetscCall(TaoComputeObjective(tao,X,f)); 246 PetscCall(TaoDefaultComputeGradient(tao,X,G,NULL)); 247 } else { 248 PetscStackPush("Tao user objective/gradient evaluation routine"); 249 PetscCall((*tao->ops->computeobjectiveandgradient)(tao,X,f,G,tao->user_objgradP)); 250 PetscStackPop; 251 } 252 PetscCall(PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL)); 253 tao->nfuncgrads++; 254 } else if (tao->ops->computeobjective && tao->ops->computegradient) { 255 PetscCall(PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL)); 256 PetscStackPush("Tao user objective evaluation routine"); 257 PetscCall((*tao->ops->computeobjective)(tao,X,f,tao->user_objP)); 258 PetscStackPop; 259 PetscCall(PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL)); 260 tao->nfuncs++; 261 PetscCall(PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL)); 262 PetscStackPush("Tao user gradient evaluation routine"); 263 PetscCall((*tao->ops->computegradient)(tao,X,G,tao->user_gradP)); 264 PetscStackPop; 265 PetscCall(PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL)); 266 tao->ngrads++; 267 } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjective() or TaoSetGradient() not set"); 268 PetscCall(PetscInfo(tao,"TAO Function evaluation: %20.19e\n",(double)(*f))); 269 PetscCall(VecLockReadPop(X)); 270 271 PetscCall(TaoTestGradient(tao,X,G)); 272 PetscFunctionReturn(0); 273 } 274 275 /*@C 276 TaoSetObjective - Sets the function evaluation routine for minimization 277 278 Logically collective on Tao 279 280 Input Parameters: 281 + tao - the Tao context 282 . func - the objective function 283 - ctx - [optional] user-defined context for private data for the function evaluation 284 routine (may be NULL) 285 286 Calling sequence of func: 287 $ func (Tao tao, Vec x, PetscReal *f, void *ctx); 288 289 + x - input vector 290 . f - function value 291 - ctx - [optional] user-defined function context 292 293 Level: beginner 294 295 .seealso: `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetObjective()` 296 @*/ 297 PetscErrorCode TaoSetObjective(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx) 298 { 299 PetscFunctionBegin; 300 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 301 if (ctx) tao->user_objP = ctx; 302 if (func) tao->ops->computeobjective = func; 303 PetscFunctionReturn(0); 304 } 305 306 /*@C 307 TaoGetObjective - Gets the function evaluation routine for minimization 308 309 Not collective 310 311 Input Parameter: 312 . tao - the Tao context 313 314 Output Parameters 315 + func - the objective function 316 - ctx - the user-defined context for private data for the function evaluation 317 318 Calling sequence of func: 319 $ func (Tao tao, Vec x, PetscReal *f, void *ctx); 320 321 + x - input vector 322 . f - function value 323 - ctx - [optional] user-defined function context 324 325 Level: beginner 326 327 .seealso: `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjective()` 328 @*/ 329 PetscErrorCode TaoGetObjective(Tao tao, PetscErrorCode (**func)(Tao, Vec, PetscReal*,void*),void **ctx) 330 { 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 333 if (func) *func = tao->ops->computeobjective; 334 if (ctx) *ctx = tao->user_objP; 335 PetscFunctionReturn(0); 336 } 337 338 /*@C 339 TaoSetResidualRoutine - Sets the residual evaluation routine for least-square applications 340 341 Logically collective on Tao 342 343 Input Parameters: 344 + tao - the Tao context 345 . func - the residual evaluation routine 346 - ctx - [optional] user-defined context for private data for the function evaluation 347 routine (may be NULL) 348 349 Calling sequence of func: 350 $ func (Tao tao, Vec x, Vec f, void *ctx); 351 352 + x - input vector 353 . f - function value vector 354 - ctx - [optional] user-defined function context 355 356 Level: beginner 357 358 .seealso: `TaoSetObjective()`, `TaoSetJacobianRoutine()` 359 @*/ 360 PetscErrorCode TaoSetResidualRoutine(Tao tao, Vec res, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx) 361 { 362 PetscFunctionBegin; 363 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 364 PetscValidHeaderSpecific(res,VEC_CLASSID,2); 365 PetscCall(PetscObjectReference((PetscObject)res)); 366 if (tao->ls_res) { 367 PetscCall(VecDestroy(&tao->ls_res)); 368 } 369 tao->ls_res = res; 370 tao->user_lsresP = ctx; 371 tao->ops->computeresidual = func; 372 373 PetscFunctionReturn(0); 374 } 375 376 /*@ 377 TaoSetResidualWeights - Give weights for the residual values. A vector can be used if only diagonal terms are used, otherwise a matrix can be give. If this function is not used, or if sigma_v and sigma_w are both NULL, then the default identity matrix will be used for weights. 378 379 Collective on Tao 380 381 Input Parameters: 382 + tao - the Tao context 383 . sigma_v - vector of weights (diagonal terms only) 384 . n - the number of weights (if using off-diagonal) 385 . rows - index list of rows for sigma_w 386 . cols - index list of columns for sigma_w 387 - vals - array of weights 388 389 Note: Either sigma_v or sigma_w (or both) should be NULL 390 391 Level: intermediate 392 393 .seealso: `TaoSetResidualRoutine()` 394 @*/ 395 PetscErrorCode TaoSetResidualWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals) 396 { 397 PetscInt i; 398 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 401 if (sigma_v) PetscValidHeaderSpecific(sigma_v,VEC_CLASSID,2); 402 PetscCall(PetscObjectReference((PetscObject)sigma_v)); 403 PetscCall(VecDestroy(&tao->res_weights_v)); 404 tao->res_weights_v = sigma_v; 405 if (vals) { 406 PetscCall(PetscFree(tao->res_weights_rows)); 407 PetscCall(PetscFree(tao->res_weights_cols)); 408 PetscCall(PetscFree(tao->res_weights_w)); 409 PetscCall(PetscMalloc1(n,&tao->res_weights_rows)); 410 PetscCall(PetscMalloc1(n,&tao->res_weights_cols)); 411 PetscCall(PetscMalloc1(n,&tao->res_weights_w)); 412 tao->res_weights_n = n; 413 for (i=0;i<n;i++) { 414 tao->res_weights_rows[i] = rows[i]; 415 tao->res_weights_cols[i] = cols[i]; 416 tao->res_weights_w[i] = vals[i]; 417 } 418 } else { 419 tao->res_weights_n = 0; 420 tao->res_weights_rows = NULL; 421 tao->res_weights_cols = NULL; 422 } 423 PetscFunctionReturn(0); 424 } 425 426 /*@ 427 TaoComputeResidual - Computes a least-squares residual vector at a given point 428 429 Collective on Tao 430 431 Input Parameters: 432 + tao - the Tao context 433 - X - input vector 434 435 Output Parameter: 436 . f - Objective vector at X 437 438 Notes: 439 TaoComputeResidual() is typically used within minimization implementations, 440 so most users would not generally call this routine themselves. 441 442 Level: advanced 443 444 .seealso: `TaoSetResidualRoutine()` 445 @*/ 446 PetscErrorCode TaoComputeResidual(Tao tao, Vec X, Vec F) 447 { 448 PetscFunctionBegin; 449 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 450 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 451 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 452 PetscCheckSameComm(tao,1,X,2); 453 PetscCheckSameComm(tao,1,F,3); 454 if (tao->ops->computeresidual) { 455 PetscCall(PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL)); 456 PetscStackPush("Tao user least-squares residual evaluation routine"); 457 PetscCall((*tao->ops->computeresidual)(tao,X,F,tao->user_lsresP)); 458 PetscStackPop; 459 PetscCall(PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL)); 460 tao->nfuncs++; 461 } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetResidualRoutine() has not been called"); 462 PetscCall(PetscInfo(tao,"TAO least-squares residual evaluation.\n")); 463 PetscFunctionReturn(0); 464 } 465 466 /*@C 467 TaoSetGradient - Sets the gradient evaluation routine for minimization 468 469 Logically collective on Tao 470 471 Input Parameters: 472 + tao - the Tao context 473 . g - [optional] the vector to internally hold the gradient computation 474 . func - the gradient function 475 - ctx - [optional] user-defined context for private data for the gradient evaluation 476 routine (may be NULL) 477 478 Calling sequence of func: 479 $ func (Tao tao, Vec x, Vec g, void *ctx); 480 481 + x - input vector 482 . g - gradient value (output) 483 - ctx - [optional] user-defined function context 484 485 Level: beginner 486 487 .seealso: `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetGradient()` 488 @*/ 489 PetscErrorCode TaoSetGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx) 490 { 491 PetscFunctionBegin; 492 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 493 if (g) { 494 PetscValidHeaderSpecific(g,VEC_CLASSID,2); 495 PetscCheckSameComm(tao,1,g,2); 496 PetscCall(PetscObjectReference((PetscObject)g)); 497 PetscCall(VecDestroy(&tao->gradient)); 498 tao->gradient = g; 499 } 500 if (func) tao->ops->computegradient = func; 501 if (ctx) tao->user_gradP = ctx; 502 PetscFunctionReturn(0); 503 } 504 505 /*@C 506 TaoGetGradient - Gets the gradient evaluation routine for minimization 507 508 Not collective 509 510 Input Parameter: 511 . tao - the Tao context 512 513 Output Parameters: 514 + g - the vector to internally hold the gradient computation 515 . func - the gradient function 516 - ctx - user-defined context for private data for the gradient evaluation routine 517 518 Calling sequence of func: 519 $ func (Tao tao, Vec x, Vec g, void *ctx); 520 521 + x - input vector 522 . g - gradient value (output) 523 - ctx - [optional] user-defined function context 524 525 Level: beginner 526 527 .seealso: `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetGradient()` 528 @*/ 529 PetscErrorCode TaoGetGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao, Vec, Vec, void*),void **ctx) 530 { 531 PetscFunctionBegin; 532 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 533 if (g) *g = tao->gradient; 534 if (func) *func = tao->ops->computegradient; 535 if (ctx) *ctx = tao->user_gradP; 536 PetscFunctionReturn(0); 537 } 538 539 /*@C 540 TaoSetObjectiveAndGradient - Sets a combined objective function and gradient evaluation routine for minimization 541 542 Logically collective on Tao 543 544 Input Parameters: 545 + tao - the Tao context 546 . g - [optional] the vector to internally hold the gradient computation 547 . func - the gradient function 548 - ctx - [optional] user-defined context for private data for the gradient evaluation 549 routine (may be NULL) 550 551 Calling sequence of func: 552 $ func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx); 553 554 + x - input vector 555 . f - objective value (output) 556 . g - gradient value (output) 557 - ctx - [optional] user-defined function context 558 559 Level: beginner 560 561 .seealso: `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetGradient()`, `TaoGetObjectiveAndGradient()` 562 @*/ 563 PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 564 { 565 PetscFunctionBegin; 566 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 567 if (g) { 568 PetscValidHeaderSpecific(g,VEC_CLASSID,2); 569 PetscCheckSameComm(tao,1,g,2); 570 PetscCall(PetscObjectReference((PetscObject)g)); 571 PetscCall(VecDestroy(&tao->gradient)); 572 tao->gradient = g; 573 } 574 if (ctx) tao->user_objgradP = ctx; 575 if (func) tao->ops->computeobjectiveandgradient = func; 576 PetscFunctionReturn(0); 577 } 578 579 /*@C 580 TaoGetObjectiveAndGradient - Gets a combined objective function and gradient evaluation routine for minimization 581 582 Not collective 583 584 Input Parameter: 585 . tao - the Tao context 586 587 Output Parameters: 588 + g - the vector to internally hold the gradient computation 589 . func - the gradient function 590 - ctx - user-defined context for private data for the gradient evaluation routine 591 592 Calling sequence of func: 593 $ func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx); 594 595 + x - input vector 596 . f - objective value (output) 597 . g - gradient value (output) 598 - ctx - [optional] user-defined function context 599 600 Level: beginner 601 602 .seealso: `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()` 603 @*/ 604 PetscErrorCode TaoGetObjectiveAndGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao, Vec, PetscReal*, Vec, void*), void **ctx) 605 { 606 PetscFunctionBegin; 607 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 608 if (g) *g = tao->gradient; 609 if (func) *func = tao->ops->computeobjectiveandgradient; 610 if (ctx) *ctx = tao->user_objgradP; 611 PetscFunctionReturn(0); 612 } 613 614 /*@ 615 TaoIsObjectiveDefined - Checks to see if the user has 616 declared an objective-only routine. Useful for determining when 617 it is appropriate to call TaoComputeObjective() or 618 TaoComputeObjectiveAndGradient() 619 620 Not collective 621 622 Input Parameter: 623 . tao - the Tao context 624 625 Output Parameter: 626 . flg - PETSC_TRUE if function routine is set by user, PETSC_FALSE otherwise 627 628 Level: developer 629 630 .seealso: `TaoSetObjective()`, `TaoIsGradientDefined()`, `TaoIsObjectiveAndGradientDefined()` 631 @*/ 632 PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg) 633 { 634 PetscFunctionBegin; 635 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 636 if (tao->ops->computeobjective == NULL) *flg = PETSC_FALSE; 637 else *flg = PETSC_TRUE; 638 PetscFunctionReturn(0); 639 } 640 641 /*@ 642 TaoIsGradientDefined - Checks to see if the user has 643 declared an objective-only routine. Useful for determining when 644 it is appropriate to call TaoComputeGradient() or 645 TaoComputeGradientAndGradient() 646 647 Not Collective 648 649 Input Parameter: 650 . tao - the Tao context 651 652 Output Parameter: 653 . flg - PETSC_TRUE if function routine is set by user, PETSC_FALSE otherwise 654 655 Level: developer 656 657 .seealso: `TaoSetGradient()`, `TaoIsObjectiveDefined()`, `TaoIsObjectiveAndGradientDefined()` 658 @*/ 659 PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg) 660 { 661 PetscFunctionBegin; 662 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 663 if (tao->ops->computegradient == NULL) *flg = PETSC_FALSE; 664 else *flg = PETSC_TRUE; 665 PetscFunctionReturn(0); 666 } 667 668 /*@ 669 TaoIsObjectiveAndGradientDefined - Checks to see if the user has 670 declared a joint objective/gradient routine. Useful for determining when 671 it is appropriate to call TaoComputeObjective() or 672 TaoComputeObjectiveAndGradient() 673 674 Not Collective 675 676 Input Parameter: 677 . tao - the Tao context 678 679 Output Parameter: 680 . flg - PETSC_TRUE if function routine is set by user, PETSC_FALSE otherwise 681 682 Level: developer 683 684 .seealso: `TaoSetObjectiveAndGradient()`, `TaoIsObjectiveDefined()`, `TaoIsGradientDefined()` 685 @*/ 686 PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg) 687 { 688 PetscFunctionBegin; 689 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 690 if (tao->ops->computeobjectiveandgradient == NULL) *flg = PETSC_FALSE; 691 else *flg = PETSC_TRUE; 692 PetscFunctionReturn(0); 693 } 694