1 #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/ 2 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 3 #include <../src/tao/bound/impls/blmvm/blmvm.h> 4 5 /*------------------------------------------------------------*/ 6 static PetscErrorCode TaoSolve_BLMVM(Tao tao) 7 { 8 PetscErrorCode ierr; 9 TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 10 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 11 PetscReal f, fold, gdx, gnorm, gnorm2; 12 PetscReal stepsize = 1.0,delta; 13 14 PetscFunctionBegin; 15 /* Project initial point onto bounds */ 16 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 17 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 18 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 19 20 21 /* Check convergence criteria */ 22 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr); 23 ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 24 25 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 26 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 27 28 tao->reason = TAO_CONTINUE_ITERATING; 29 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 30 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);CHKERRQ(ierr); 31 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 32 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 33 34 /* Set counter for gradient/reset steps */ 35 if (!blmP->recycle) { 36 blmP->grad = 0; 37 blmP->reset = 0; 38 ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr); 39 } 40 41 /* Have not converged; continue with Newton method */ 42 while (tao->reason == TAO_CONTINUE_ITERATING) { 43 /* Call general purpose update function */ 44 if (tao->ops->update) { 45 ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 46 } 47 /* Compute direction */ 48 gnorm2 = gnorm*gnorm; 49 if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON; 50 if (f == 0.0) { 51 delta = 2.0 / gnorm2; 52 } else { 53 delta = 2.0 * PetscAbsScalar(f) / gnorm2; 54 } 55 ierr = MatSymBrdnSetDelta(blmP->M, delta);CHKERRQ(ierr); 56 ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 57 ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 58 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 59 60 /* Check for success (descent direction) */ 61 ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr); 62 if (gdx <= 0) { 63 /* Step is not descent or solve was not successful 64 Use steepest descent direction (scaled) */ 65 ++blmP->grad; 66 67 ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr); 68 ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 69 ierr = MatSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 70 } 71 ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr); 72 73 /* Perform the linesearch */ 74 fold = f; 75 ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr); 76 ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr); 77 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 78 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr); 79 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 80 81 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 82 /* Linesearch failed 83 Reset factors and use scaled (projected) gradient step */ 84 ++blmP->reset; 85 86 f = fold; 87 ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr); 88 ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr); 89 90 ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr); 91 ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 92 ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 93 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 94 95 /* This may be incorrect; linesearch has values for stepmax and stepmin 96 that should be reset. */ 97 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 98 ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr); 99 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 100 101 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 102 tao->reason = TAO_DIVERGED_LS_FAILURE; 103 break; 104 } 105 } 106 107 /* Check for converged */ 108 ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 109 ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 110 111 112 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number"); 113 tao->niter++; 114 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 115 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);CHKERRQ(ierr); 116 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 117 } 118 PetscFunctionReturn(0); 119 } 120 121 static PetscErrorCode TaoSetup_BLMVM(Tao tao) 122 { 123 TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 124 PetscErrorCode ierr; 125 126 PetscFunctionBegin; 127 /* Existence of tao->solution checked in TaoSetup() */ 128 ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr); 129 ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr); 130 ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr); 131 132 if (!tao->stepdirection) { 133 ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 134 } 135 if (!tao->gradient) { 136 ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 137 } 138 if (!tao->XL) { 139 ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr); 140 ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr); 141 } 142 if (!tao->XU) { 143 ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr); 144 ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr); 145 } 146 /* Allocate matrix for the limited memory approximation */ 147 ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr); 148 149 /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 150 if (blmP->H0) { 151 ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr); 152 } 153 PetscFunctionReturn(0); 154 } 155 156 /* ---------------------------------------------------------- */ 157 static PetscErrorCode TaoDestroy_BLMVM(Tao tao) 158 { 159 TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 160 PetscErrorCode ierr; 161 162 PetscFunctionBegin; 163 if (tao->setupcalled) { 164 ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr); 165 ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr); 166 ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr); 167 } 168 ierr = MatDestroy(&blmP->M);CHKERRQ(ierr); 169 if (blmP->H0) { 170 PetscObjectDereference((PetscObject)blmP->H0); 171 } 172 ierr = PetscFree(tao->data);CHKERRQ(ierr); 173 PetscFunctionReturn(0); 174 } 175 176 /*------------------------------------------------------------*/ 177 static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao) 178 { 179 TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 180 PetscErrorCode ierr; 181 PetscBool is_spd; 182 183 PetscFunctionBegin; 184 ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr); 185 ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr); 186 ierr = PetscOptionsTail();CHKERRQ(ierr); 187 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 188 ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr); 189 ierr = MatGetOption(blmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr); 190 if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 191 PetscFunctionReturn(0); 192 } 193 194 195 /*------------------------------------------------------------*/ 196 static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) 197 { 198 TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 199 PetscBool isascii; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 204 if (isascii) { 205 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr); 206 ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 207 ierr = MatView(lmP->M, viewer);CHKERRQ(ierr); 208 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 209 } 210 PetscFunctionReturn(0); 211 } 212 213 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 214 { 215 TAO_BLMVM *blm = (TAO_BLMVM *) tao->data; 216 PetscErrorCode ierr; 217 218 PetscFunctionBegin; 219 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 220 PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 221 PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 222 if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 223 224 ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr); 225 ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr); 226 ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 227 ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 228 229 ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr); 230 ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr); 231 ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr); 232 PetscFunctionReturn(0); 233 } 234 235 /* ---------------------------------------------------------- */ 236 /*MC 237 TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 238 for nonlinear minimization with bound constraints. It is an extension 239 of TAOLMVM 240 241 Options Database Keys: 242 . -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls 243 244 Level: beginner 245 M*/ 246 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 247 { 248 TAO_BLMVM *blmP; 249 const char *morethuente_type = TAOLINESEARCHMT; 250 PetscErrorCode ierr; 251 252 PetscFunctionBegin; 253 tao->ops->setup = TaoSetup_BLMVM; 254 tao->ops->solve = TaoSolve_BLMVM; 255 tao->ops->view = TaoView_BLMVM; 256 tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 257 tao->ops->destroy = TaoDestroy_BLMVM; 258 tao->ops->computedual = TaoComputeDual_BLMVM; 259 260 ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr); 261 blmP->H0 = NULL; 262 blmP->recycle = PETSC_FALSE; 263 tao->data = (void*)blmP; 264 265 /* Override default settings (unless already changed) */ 266 if (!tao->max_it_changed) tao->max_it = 2000; 267 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 268 269 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 270 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 271 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 272 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 273 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 274 275 ierr = KSPInitializePackage();CHKERRQ(ierr); 276 ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr); 277 ierr = MatSetType(blmP->M, MATLMVMBFGS);CHKERRQ(ierr); 278 ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr); 279 ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr); 280 PetscFunctionReturn(0); 281 } 282 283 /*@ 284 TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls. 285 286 Input Parameters: 287 + tao - the Tao solver context 288 - flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE) 289 290 Level: intermediate 291 @*/ 292 PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg) 293 { 294 TAO_LMVM *lmP; 295 TAO_BLMVM *blmP; 296 TaoType type; 297 PetscBool is_lmvm, is_blmvm; 298 PetscErrorCode ierr; 299 300 PetscFunctionBegin; 301 ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 302 ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 303 ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 304 305 if (is_lmvm) { 306 lmP = (TAO_LMVM *)tao->data; 307 lmP->recycle = flg; 308 } else if (is_blmvm) { 309 blmP = (TAO_BLMVM *)tao->data; 310 blmP->recycle = flg; 311 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 312 PetscFunctionReturn(0); 313 } 314 315 /*@ 316 TaoLMVMSetH0 - Set the initial Hessian for the QN approximation 317 318 Input Parameters: 319 + tao - the Tao solver context 320 - H0 - Mat object for the initial Hessian 321 322 Level: advanced 323 324 .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP() 325 @*/ 326 PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) 327 { 328 TAO_LMVM *lmP; 329 TAO_BLMVM *blmP; 330 TaoType type; 331 PetscBool is_lmvm, is_blmvm; 332 PetscErrorCode ierr; 333 334 PetscFunctionBegin; 335 ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 336 ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 337 ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 338 339 if (is_lmvm) { 340 lmP = (TAO_LMVM *)tao->data; 341 ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 342 lmP->H0 = H0; 343 } else if (is_blmvm) { 344 blmP = (TAO_BLMVM *)tao->data; 345 ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 346 blmP->H0 = H0; 347 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 348 PetscFunctionReturn(0); 349 } 350 351 /*@ 352 TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian 353 354 Input Parameters: 355 . tao - the Tao solver context 356 357 Output Parameters: 358 . H0 - Mat object for the initial Hessian 359 360 Level: advanced 361 362 .seealso: TaoLMVMSetH0(), TaoLMVMGetH0KSP() 363 @*/ 364 PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) 365 { 366 TAO_LMVM *lmP; 367 TAO_BLMVM *blmP; 368 TaoType type; 369 PetscBool is_lmvm, is_blmvm; 370 Mat M; 371 372 PetscErrorCode ierr; 373 374 PetscFunctionBegin; 375 ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 376 ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 377 ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 378 379 if (is_lmvm) { 380 lmP = (TAO_LMVM *)tao->data; 381 M = lmP->M; 382 } else if (is_blmvm) { 383 blmP = (TAO_BLMVM *)tao->data; 384 M = blmP->M; 385 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 386 ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr); 387 PetscFunctionReturn(0); 388 } 389 390 /*@ 391 TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian 392 393 Input Parameters: 394 . tao - the Tao solver context 395 396 Output Parameters: 397 . ksp - KSP solver context for the initial Hessian 398 399 Level: advanced 400 401 .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP() 402 @*/ 403 PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) 404 { 405 TAO_LMVM *lmP; 406 TAO_BLMVM *blmP; 407 TaoType type; 408 PetscBool is_lmvm, is_blmvm; 409 Mat M; 410 PetscErrorCode ierr; 411 412 ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 413 ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 414 ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 415 416 if (is_lmvm) { 417 lmP = (TAO_LMVM *)tao->data; 418 M = lmP->M; 419 } else if (is_blmvm) { 420 blmP = (TAO_BLMVM *)tao->data; 421 M = blmP->M; 422 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 423 ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr); 424 PetscFunctionReturn(0); 425 } 426