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