1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "TaoSetInitialVector" 5 /*@ 6 TaoSetInitialVector - Sets the initial guess for the solve 7 8 Logically collective on Tao 9 10 Input Parameters: 11 + tao - the Tao context 12 - x0 - the initial guess 13 14 Level: beginner 15 .seealso: TaoCreate(), TaoSolve() 16 @*/ 17 18 PetscErrorCode TaoSetInitialVector(Tao tao, Vec x0) 19 { 20 PetscErrorCode ierr; 21 22 PetscFunctionBegin; 23 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 24 if (x0) { 25 PetscValidHeaderSpecific(x0,VEC_CLASSID,2); 26 PetscObjectReference((PetscObject)x0); 27 } 28 ierr = VecDestroy(&tao->solution);CHKERRQ(ierr); 29 tao->solution = x0; 30 PetscFunctionReturn(0); 31 } 32 33 #undef __FUNCT__ 34 #define __FUNCT__ "TaoComputeGradient" 35 /*@ 36 TaoComputeGradient - Computes the gradient of the objective function 37 38 Collective on Tao 39 40 Input Parameters: 41 + tao - the Tao context 42 - X - input vector 43 44 Output Parameter: 45 . G - gradient vector 46 47 Notes: TaoComputeGradient() is typically used within minimization implementations, 48 so most users would not generally call this routine themselves. 49 50 Level: advanced 51 52 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetGradientRoutine() 53 @*/ 54 PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G) 55 { 56 PetscErrorCode ierr; 57 PetscReal dummy; 58 59 PetscFunctionBegin; 60 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 61 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 62 PetscValidHeaderSpecific(G,VEC_CLASSID,2); 63 PetscCheckSameComm(tao,1,X,2); 64 PetscCheckSameComm(tao,1,G,3); 65 if (tao->ops->computegradient) { 66 ierr = PetscLogEventBegin(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr); 67 PetscStackPush("Tao user gradient evaluation routine"); 68 ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr); 69 PetscStackPop; 70 ierr = PetscLogEventEnd(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr); 71 tao->ngrads++; 72 } else if (tao->ops->computeobjectiveandgradient) { 73 ierr = PetscLogEventBegin(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr); 74 PetscStackPush("Tao user objective/gradient evaluation routine"); 75 ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,&dummy,G,tao->user_objgradP);CHKERRQ(ierr); 76 PetscStackPop; 77 ierr = PetscLogEventEnd(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr); 78 tao->nfuncgrads++; 79 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetGradientRoutine() has not been called"); 80 PetscFunctionReturn(0); 81 } 82 83 #undef __FUNCT__ 84 #define __FUNCT__ "TaoComputeObjective" 85 /*@ 86 TaoComputeObjective - Computes the objective function value at a given point 87 88 Collective on Tao 89 90 Input Parameters: 91 + tao - the Tao context 92 - X - input vector 93 94 Output Parameter: 95 . f - Objective value at X 96 97 Notes: TaoComputeObjective() is typically used within minimization implementations, 98 so most users would not generally call this routine themselves. 99 100 Level: advanced 101 102 .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine() 103 @*/ 104 PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f) 105 { 106 PetscErrorCode ierr; 107 Vec temp; 108 109 PetscFunctionBegin; 110 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 111 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 112 PetscCheckSameComm(tao,1,X,2); 113 if (tao->ops->computeobjective) { 114 ierr = PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr); 115 PetscStackPush("Tao user objective evaluation routine"); 116 ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr); 117 PetscStackPop; 118 ierr = PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr); 119 tao->nfuncs++; 120 } else if (tao->ops->computeobjectiveandgradient) { 121 ierr = PetscInfo(tao,"Duplicating variable vector in order to call func/grad routine\n");CHKERRQ(ierr); 122 ierr = VecDuplicate(X,&temp);CHKERRQ(ierr); 123 ierr = PetscLogEventBegin(Tao_ObjGradientEval,tao,X,NULL,NULL);CHKERRQ(ierr); 124 PetscStackPush("Tao user objective/gradient evaluation routine"); 125 ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,temp,tao->user_objgradP);CHKERRQ(ierr); 126 PetscStackPop; 127 ierr = PetscLogEventEnd(Tao_ObjGradientEval,tao,X,NULL,NULL);CHKERRQ(ierr); 128 ierr = VecDestroy(&temp);CHKERRQ(ierr); 129 tao->nfuncgrads++; 130 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() has not been called"); 131 ierr = PetscInfo1(tao,"TAO Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 132 PetscFunctionReturn(0); 133 } 134 135 #undef __FUNCT__ 136 #define __FUNCT__ "TaoComputeObjectiveAndGradient" 137 /*@ 138 TaoComputeObjectiveAndGradient - Computes the objective function value at a given point 139 140 Collective on Tao 141 142 Input Parameters: 143 + tao - the Tao context 144 - X - input vector 145 146 Output Parameter: 147 + f - Objective value at X 148 - g - Gradient vector at X 149 150 Notes: TaoComputeObjectiveAndGradient() is typically used within minimization implementations, 151 so most users would not generally call this routine themselves. 152 153 Level: advanced 154 155 .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine() 156 @*/ 157 PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G) 158 { 159 PetscErrorCode ierr; 160 161 PetscFunctionBegin; 162 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 163 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 164 PetscValidHeaderSpecific(G,VEC_CLASSID,4); 165 PetscCheckSameComm(tao,1,X,2); 166 PetscCheckSameComm(tao,1,G,4); 167 if (tao->ops->computeobjectiveandgradient) { 168 ierr = PetscLogEventBegin(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr); 169 PetscStackPush("Tao user objective/gradient evaluation routine"); 170 ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,G,tao->user_objgradP);CHKERRQ(ierr); 171 PetscStackPop; 172 if (tao->ops->computegradient == TaoDefaultComputeGradient) { 173 /* Overwrite gradient with finite difference gradient */ 174 ierr = TaoDefaultComputeGradient(tao,X,G,tao->user_objgradP);CHKERRQ(ierr); 175 } 176 ierr = PetscLogEventEnd(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr); 177 tao->nfuncgrads++; 178 } else if (tao->ops->computeobjective && tao->ops->computegradient) { 179 ierr = PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr); 180 PetscStackPush("Tao user objective evaluation routine"); 181 ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr); 182 PetscStackPop; 183 ierr = PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr); 184 tao->nfuncs++; 185 ierr = PetscLogEventBegin(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr); 186 PetscStackPush("Tao user gradient evaluation routine"); 187 ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr); 188 PetscStackPop; 189 ierr = PetscLogEventEnd(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr); 190 tao->ngrads++; 191 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() or TaoSetGradientRoutine() not set"); 192 ierr = PetscInfo1(tao,"TAO Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 193 PetscFunctionReturn(0); 194 } 195 196 #undef __FUNCT__ 197 #define __FUNCT__ "TaoSetObjectiveRoutine" 198 /*@C 199 TaoSetObjectiveRoutine - Sets the function evaluation routine for minimization 200 201 Logically collective on Tao 202 203 Input Parameter: 204 + tao - the Tao context 205 . func - the objective function 206 - ctx - [optional] user-defined context for private data for the function evaluation 207 routine (may be NULL) 208 209 Calling sequence of func: 210 $ func (Tao tao, Vec x, PetscReal *f, void *ctx); 211 212 + x - input vector 213 . f - function value 214 - ctx - [optional] user-defined function context 215 216 Level: beginner 217 218 .seealso: TaoSetGradientRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine() 219 @*/ 220 PetscErrorCode TaoSetObjectiveRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx) 221 { 222 PetscFunctionBegin; 223 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 224 tao->user_objP = ctx; 225 tao->ops->computeobjective = func; 226 PetscFunctionReturn(0); 227 } 228 229 #undef __FUNCT__ 230 #define __FUNCT__ "TaoSetSeparableObjectiveRoutine" 231 /*@C 232 TaoSetSeparableObjectiveRoutine - Sets the function evaluation routine for least-square applications 233 234 Logically collective on Tao 235 236 Input Parameter: 237 + tao - the Tao context 238 . func - the objective function evaluation routine 239 - ctx - [optional] user-defined context for private data for the function evaluation 240 routine (may be NULL) 241 242 Calling sequence of func: 243 $ func (Tao tao, Vec x, Vec f, void *ctx); 244 245 + x - input vector 246 . f - function value vector 247 - ctx - [optional] user-defined function context 248 249 Level: beginner 250 251 .seealso: TaoSetObjectiveRoutine(), TaoSetJacobianRoutine() 252 @*/ 253 PetscErrorCode TaoSetSeparableObjectiveRoutine(Tao tao, Vec sepobj, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx) 254 { 255 PetscFunctionBegin; 256 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 257 PetscValidHeaderSpecific(sepobj, VEC_CLASSID,2); 258 tao->user_sepobjP = ctx; 259 tao->sep_objective = sepobj; 260 tao->ops->computeseparableobjective = func; 261 PetscFunctionReturn(0); 262 } 263 264 #undef __FUNCT__ 265 #define __FUNCT__ "TaoSetSeparableObjectiveWeights" 266 /*@ 267 TaoSetSeparableObjectiveWeights - Give weights for the separable objective 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. 268 269 Collective on Tao 270 271 Input Parameters: 272 + tao - the Tao context 273 . sigma_v - vector of weights (diagonal terms only) 274 . n - the number of weights (if using off-diagonal) 275 . rows - index list of rows for sigma_w 276 . cols - index list of columns for sigma_w 277 - vals - array of weights 278 279 280 281 Note: Either sigma_v or sigma_w (or both) should be NULL 282 283 Level: intermediate 284 285 .seealso: TaoSetSeparableObjectiveRoutine() 286 @*/ 287 PetscErrorCode TaoSetSeparableObjectiveWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals) 288 { 289 PetscErrorCode ierr; 290 PetscInt i; 291 PetscFunctionBegin; 292 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 293 ierr = VecDestroy(&tao->sep_weights_v);CHKERRQ(ierr); 294 tao->sep_weights_v=sigma_v; 295 if (sigma_v) { 296 ierr = PetscObjectReference((PetscObject)sigma_v);CHKERRQ(ierr); 297 } 298 if (vals) { 299 if (tao->sep_weights_n) { 300 ierr = PetscFree(tao->sep_weights_rows);CHKERRQ(ierr); 301 ierr = PetscFree(tao->sep_weights_cols);CHKERRQ(ierr); 302 ierr = PetscFree(tao->sep_weights_w);CHKERRQ(ierr); 303 } 304 ierr = PetscMalloc1(n,&tao->sep_weights_rows);CHKERRQ(ierr); 305 ierr = PetscMalloc1(n,&tao->sep_weights_cols);CHKERRQ(ierr); 306 ierr = PetscMalloc1(n,&tao->sep_weights_w);CHKERRQ(ierr); 307 tao->sep_weights_n=n; 308 for (i=0;i<n;i++) { 309 tao->sep_weights_rows[i]=rows[i]; 310 tao->sep_weights_cols[i]=cols[i]; 311 tao->sep_weights_w[i]=vals[i]; 312 } 313 } else { 314 tao->sep_weights_n=0; 315 tao->sep_weights_rows=0; 316 tao->sep_weights_cols=0; 317 } 318 PetscFunctionReturn(0); 319 } 320 #undef __FUNCT__ 321 #define __FUNCT__ "TaoComputeSeparableObjective" 322 /*@ 323 TaoComputeSeparableObjective - Computes a separable objective function vector at a given point (for least-square applications) 324 325 Collective on Tao 326 327 Input Parameters: 328 + tao - the Tao context 329 - X - input vector 330 331 Output Parameter: 332 . f - Objective vector at X 333 334 Notes: TaoComputeSeparableObjective() is typically used within minimization implementations, 335 so most users would not generally call this routine themselves. 336 337 Level: advanced 338 339 .seealso: TaoSetSeparableObjectiveRoutine() 340 @*/ 341 PetscErrorCode TaoComputeSeparableObjective(Tao tao, Vec X, Vec F) 342 { 343 PetscErrorCode ierr; 344 345 PetscFunctionBegin; 346 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 347 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 348 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 349 PetscCheckSameComm(tao,1,X,2); 350 PetscCheckSameComm(tao,1,F,3); 351 if (tao->ops->computeseparableobjective) { 352 ierr = PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr); 353 PetscStackPush("Tao user separable objective evaluation routine"); 354 ierr = (*tao->ops->computeseparableobjective)(tao,X,F,tao->user_sepobjP);CHKERRQ(ierr); 355 PetscStackPop; 356 ierr = PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr); 357 tao->nfuncs++; 358 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetSeparableObjectiveRoutine() has not been called"); 359 ierr = PetscInfo(tao,"TAO separable function evaluation.\n");CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "TaoSetGradientRoutine" 365 /*@C 366 TaoSetGradientRoutine - Sets the gradient evaluation routine for minimization 367 368 Logically collective on Tao 369 370 Input Parameter: 371 + tao - the Tao context 372 . func - the gradient function 373 - ctx - [optional] user-defined context for private data for the gradient evaluation 374 routine (may be NULL) 375 376 Calling sequence of func: 377 $ func (Tao tao, Vec x, Vec g, void *ctx); 378 379 + x - input vector 380 . g - gradient value (output) 381 - ctx - [optional] user-defined function context 382 383 Level: beginner 384 385 .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine() 386 @*/ 387 PetscErrorCode TaoSetGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx) 388 { 389 PetscFunctionBegin; 390 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 391 tao->user_gradP = ctx; 392 tao->ops->computegradient = func; 393 PetscFunctionReturn(0); 394 } 395 396 397 #undef __FUNCT__ 398 #define __FUNCT__ "TaoSetObjectiveAndGradientRoutine" 399 /*@C 400 TaoSetObjectiveAndGradientRoutine - Sets a combined objective function and gradient evaluation routine for minimization 401 402 Logically collective on Tao 403 404 Input Parameter: 405 + tao - the Tao context 406 . func - the gradient function 407 - ctx - [optional] user-defined context for private data for the gradient evaluation 408 routine (may be NULL) 409 410 Calling sequence of func: 411 $ func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx); 412 413 + x - input vector 414 . f - objective value (output) 415 . g - gradient value (output) 416 - ctx - [optional] user-defined function context 417 418 Level: beginner 419 420 .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine() 421 @*/ 422 PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void*), void *ctx) 423 { 424 PetscFunctionBegin; 425 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 426 tao->user_objgradP = ctx; 427 tao->ops->computeobjectiveandgradient = func; 428 PetscFunctionReturn(0); 429 } 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "TaoIsObjectiveDefined" 433 /*@ 434 TaoIsObjectiveDefined -- Checks to see if the user has 435 declared an objective-only routine. Useful for determining when 436 it is appropriate to call TaoComputeObjective() or 437 TaoComputeObjectiveAndGradient() 438 439 Collective on Tao 440 441 Input Parameter: 442 + tao - the Tao context 443 - ctx - PETSC_TRUE if objective function routine is set by user, 444 PETSC_FALSE otherwise 445 Level: developer 446 447 .seealso: TaoSetObjectiveRoutine(), TaoIsGradientDefined(), TaoIsObjectiveAndGradientDefined() 448 @*/ 449 PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg) 450 { 451 PetscFunctionBegin; 452 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 453 if (tao->ops->computeobjective == 0) *flg = PETSC_FALSE; 454 else *flg = PETSC_TRUE; 455 PetscFunctionReturn(0); 456 } 457 458 #undef __FUNCT__ 459 #define __FUNCT__ "TaoIsGradientDefined" 460 /*@ 461 TaoIsGradientDefined -- Checks to see if the user has 462 declared an objective-only routine. Useful for determining when 463 it is appropriate to call TaoComputeGradient() or 464 TaoComputeGradientAndGradient() 465 466 Not Collective 467 468 Input Parameter: 469 + tao - the Tao context 470 - ctx - PETSC_TRUE if gradient routine is set by user, PETSC_FALSE otherwise 471 Level: developer 472 473 .seealso: TaoSetGradientRoutine(), TaoIsObjectiveDefined(), TaoIsObjectiveAndGradientDefined() 474 @*/ 475 PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg) 476 { 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 479 if (tao->ops->computegradient == 0) *flg = PETSC_FALSE; 480 else *flg = PETSC_TRUE; 481 PetscFunctionReturn(0); 482 } 483 484 485 #undef __FUNCT__ 486 #define __FUNCT__ "TaoIsObjectiveAndGradientDefined" 487 /*@ 488 TaoIsObjectiveAndGradientDefined -- Checks to see if the user has 489 declared a joint objective/gradient routine. Useful for determining when 490 it is appropriate to call TaoComputeObjective() or 491 TaoComputeObjectiveAndGradient() 492 493 Not Collective 494 495 Input Parameter: 496 + tao - the Tao context 497 - ctx - PETSC_TRUE if objective/gradient routine is set by user, PETSC_FALSE otherwise 498 Level: developer 499 500 .seealso: TaoSetObjectiveAndGradientRoutine(), TaoIsObjectiveDefined(), TaoIsGradientDefined() 501 @*/ 502 PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg) 503 { 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 506 if (tao->ops->computeobjectiveandgradient == 0) *flg = PETSC_FALSE; 507 else *flg = PETSC_TRUE; 508 PetscFunctionReturn(0); 509 } 510 511 512 513