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