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 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(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 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(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 = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 185 ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr); 186 ierr = MatGetOption(blmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr); 187 if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 188 PetscFunctionReturn(0); 189 } 190 191 /*------------------------------------------------------------*/ 192 static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) 193 { 194 TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 195 PetscBool isascii; 196 PetscErrorCode ierr; 197 198 PetscFunctionBegin; 199 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 200 if (isascii) { 201 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr); 202 ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 203 ierr = MatView(lmP->M, viewer);CHKERRQ(ierr); 204 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 205 } 206 PetscFunctionReturn(0); 207 } 208 209 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 210 { 211 TAO_BLMVM *blm = (TAO_BLMVM *) tao->data; 212 PetscErrorCode ierr; 213 214 PetscFunctionBegin; 215 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 216 PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 217 PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 218 if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 219 220 ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr); 221 ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr); 222 ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 223 ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 224 225 ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr); 226 ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr); 227 ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 /* ---------------------------------------------------------- */ 232 /*MC 233 TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 234 for nonlinear minimization with bound constraints. It is an extension 235 of TAOLMVM 236 237 Options Database Keys: 238 . -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls 239 240 Level: beginner 241 M*/ 242 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 243 { 244 TAO_BLMVM *blmP; 245 const char *morethuente_type = TAOLINESEARCHMT; 246 PetscErrorCode ierr; 247 248 PetscFunctionBegin; 249 tao->ops->setup = TaoSetup_BLMVM; 250 tao->ops->solve = TaoSolve_BLMVM; 251 tao->ops->view = TaoView_BLMVM; 252 tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 253 tao->ops->destroy = TaoDestroy_BLMVM; 254 tao->ops->computedual = TaoComputeDual_BLMVM; 255 256 ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr); 257 blmP->H0 = NULL; 258 blmP->recycle = PETSC_FALSE; 259 tao->data = (void*)blmP; 260 261 /* Override default settings (unless already changed) */ 262 if (!tao->max_it_changed) tao->max_it = 2000; 263 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 264 265 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 266 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 267 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 268 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 269 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);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 ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278 279 /*@ 280 TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls. 281 282 Input Parameters: 283 + tao - the Tao solver context 284 - flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE) 285 286 Level: intermediate 287 @*/ 288 PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg) 289 { 290 TAO_LMVM *lmP; 291 TAO_BLMVM *blmP; 292 PetscBool is_lmvm, is_blmvm; 293 PetscErrorCode ierr; 294 295 PetscFunctionBegin; 296 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);CHKERRQ(ierr); 297 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);CHKERRQ(ierr); 298 if (is_lmvm) { 299 lmP = (TAO_LMVM *)tao->data; 300 lmP->recycle = flg; 301 } else if (is_blmvm) { 302 blmP = (TAO_BLMVM *)tao->data; 303 blmP->recycle = flg; 304 } 305 PetscFunctionReturn(0); 306 } 307 308 /*@ 309 TaoLMVMSetH0 - Set the initial Hessian for the QN approximation 310 311 Input Parameters: 312 + tao - the Tao solver context 313 - H0 - Mat object for the initial Hessian 314 315 Level: advanced 316 317 .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP() 318 @*/ 319 PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) 320 { 321 TAO_LMVM *lmP; 322 TAO_BLMVM *blmP; 323 PetscBool is_lmvm, is_blmvm; 324 PetscErrorCode ierr; 325 326 PetscFunctionBegin; 327 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);CHKERRQ(ierr); 328 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);CHKERRQ(ierr); 329 if (is_lmvm) { 330 lmP = (TAO_LMVM *)tao->data; 331 ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 332 lmP->H0 = H0; 333 } else if (is_blmvm) { 334 blmP = (TAO_BLMVM *)tao->data; 335 ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 336 blmP->H0 = H0; 337 } 338 PetscFunctionReturn(0); 339 } 340 341 /*@ 342 TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian 343 344 Input Parameters: 345 . tao - the Tao solver context 346 347 Output Parameters: 348 . H0 - Mat object for the initial Hessian 349 350 Level: advanced 351 352 .seealso: TaoLMVMSetH0(), TaoLMVMGetH0KSP() 353 @*/ 354 PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) 355 { 356 TAO_LMVM *lmP; 357 TAO_BLMVM *blmP; 358 PetscBool is_lmvm, is_blmvm; 359 Mat M; 360 PetscErrorCode ierr; 361 362 PetscFunctionBegin; 363 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);CHKERRQ(ierr); 364 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);CHKERRQ(ierr); 365 if (is_lmvm) { 366 lmP = (TAO_LMVM *)tao->data; 367 M = lmP->M; 368 } else if (is_blmvm) { 369 blmP = (TAO_BLMVM *)tao->data; 370 M = blmP->M; 371 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 372 ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr); 373 PetscFunctionReturn(0); 374 } 375 376 /*@ 377 TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian 378 379 Input Parameters: 380 . tao - the Tao solver context 381 382 Output Parameters: 383 . ksp - KSP solver context for the initial Hessian 384 385 Level: advanced 386 387 .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP() 388 @*/ 389 PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) 390 { 391 TAO_LMVM *lmP; 392 TAO_BLMVM *blmP; 393 PetscBool is_lmvm, is_blmvm; 394 Mat M; 395 PetscErrorCode ierr; 396 397 PetscFunctionBegin; 398 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);CHKERRQ(ierr); 399 ierr = PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);CHKERRQ(ierr); 400 if (is_lmvm) { 401 lmP = (TAO_LMVM *)tao->data; 402 M = lmP->M; 403 } else if (is_blmvm) { 404 blmP = (TAO_BLMVM *)tao->data; 405 M = blmP->M; 406 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 407 ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr); 408 PetscFunctionReturn(0); 409 } 410