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