1 #include <petsctaolinesearch.h> 2 #include <../src/tao/matrix/lmvmmat.h> 3 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 4 5 #define LMVM_BFGS 0 6 #define LMVM_SCALED_GRADIENT 1 7 #define LMVM_GRADIENT 2 8 9 static PetscErrorCode TaoSolve_LMVM(Tao tao) 10 { 11 TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 12 PetscReal f, fold, gdx, gnorm; 13 PetscReal step = 1.0; 14 PetscReal delta; 15 PetscErrorCode ierr; 16 PetscInt stepType, nupdates; 17 PetscBool recycle; 18 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 19 20 PetscFunctionBegin; 21 22 if (tao->XL || tao->XU || tao->ops->computebounds) { 23 ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr); 24 } 25 26 /* Check convergence criteria */ 27 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 28 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 29 30 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 31 32 tao->reason = TAO_CONTINUE_ITERATING; 33 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 34 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 35 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 36 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 37 38 /* Set initial scaling for the function */ 39 if (f != 0.0) { 40 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 41 } else { 42 delta = 2.0 / (gnorm*gnorm); 43 } 44 ierr = MatLMVMSetDelta(lmP->M,delta);CHKERRQ(ierr); 45 46 /* Set counter for gradient/reset steps */ 47 ierr = MatLMVMGetRecycleFlag(lmP->M, &recycle);CHKERRQ(ierr); 48 if (!recycle) { 49 lmP->bfgs = 0; 50 lmP->sgrad = 0; 51 lmP->grad = 0; 52 ierr = MatLMVMReset(lmP->M); CHKERRQ(ierr); 53 } 54 55 /* Have not converged; continue with Newton method */ 56 while (tao->reason == TAO_CONTINUE_ITERATING) { 57 /* Compute direction */ 58 ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 59 ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 60 61 /* Check for success (descent direction) */ 62 ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr); 63 if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 64 /* Step is not descent or direction produced not a number 65 We can assert bfgsUpdates > 1 in this case because 66 the first solve produces the scaled gradient direction, 67 which is guaranteed to be descent 68 69 Use steepest descent direction (scaled) 70 */ 71 72 if (f != 0.0) { 73 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 74 } else { 75 delta = 2.0 / (gnorm*gnorm); 76 } 77 ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 78 ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 79 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 80 ierr = MatLMVMSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr); 81 82 /* On a reset, the direction cannot be not a number; it is a 83 scaled gradient step. No need to check for this condition. */ 84 ++lmP->sgrad; 85 stepType = LMVM_SCALED_GRADIENT; 86 } else { 87 ierr = MatLMVMGetUpdates(lmP->M, &nupdates); CHKERRQ(ierr); 88 if (1 == nupdates) { 89 /* The first BFGS direction is always the scaled gradient */ 90 ++lmP->sgrad; 91 stepType = LMVM_SCALED_GRADIENT; 92 } else { 93 ++lmP->bfgs; 94 stepType = LMVM_BFGS; 95 } 96 } 97 ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 98 99 /* Perform the linesearch */ 100 fold = f; 101 ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 102 ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 103 104 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr); 105 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 106 107 while (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_GRADIENT)) { 108 /* Reset factors and use scaled gradient step */ 109 f = fold; 110 ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 111 ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 112 113 switch(stepType) { 114 case LMVM_BFGS: 115 /* Failed to obtain acceptable iterate with BFGS step */ 116 /* Attempt to use the scaled gradient direction */ 117 118 if (f != 0.0) { 119 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 120 } else { 121 delta = 2.0 / (gnorm*gnorm); 122 } 123 ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 124 ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 125 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 126 ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 127 128 /* On a reset, the direction cannot be not a number; it is a 129 scaled gradient step. No need to check for this condition. */ 130 --lmP->bfgs; 131 ++lmP->sgrad; 132 stepType = LMVM_SCALED_GRADIENT; 133 break; 134 135 case LMVM_SCALED_GRADIENT: 136 /* The scaled gradient step did not produce a new iterate; 137 attempt to use the gradient direction. 138 Need to make sure we are not using a different diagonal scaling */ 139 ierr = MatLMVMSetDelta(lmP->M, 1.0);CHKERRQ(ierr); 140 ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 141 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 142 ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 143 144 --lmP->sgrad; 145 ++lmP->grad; 146 stepType = LMVM_GRADIENT; 147 break; 148 } 149 ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 150 151 /* Perform the linesearch */ 152 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr); 153 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 154 } 155 156 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 157 /* Failed to find an improving point */ 158 f = fold; 159 ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 160 ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 161 step = 0.0; 162 tao->reason = TAO_DIVERGED_LS_FAILURE; 163 } 164 165 /* Check for termination */ 166 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 167 168 tao->niter++; 169 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 170 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 171 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 172 } 173 PetscFunctionReturn(0); 174 } 175 176 static PetscErrorCode TaoSetUp_LMVM(Tao tao) 177 { 178 TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 179 PetscInt n,N; 180 PetscErrorCode ierr; 181 KSP H0ksp; 182 183 PetscFunctionBegin; 184 /* Existence of tao->solution checked in TaoSetUp() */ 185 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 186 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 187 if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 188 if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 189 if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 190 191 /* Create matrix for the limited memory approximation */ 192 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 193 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 194 ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&lmP->M);CHKERRQ(ierr); 195 ierr = MatLMVMAllocateVectors(lmP->M,tao->solution);CHKERRQ(ierr); 196 197 /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 198 if (lmP->H0) { 199 const char *prefix; 200 PC H0pc; 201 202 ierr = MatLMVMSetH0(lmP->M, lmP->H0);CHKERRQ(ierr); 203 ierr = MatLMVMGetH0KSP(lmP->M, &H0ksp);CHKERRQ(ierr); 204 205 ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr); 206 ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr); 207 ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0ksp, "tao_h0_");CHKERRQ(ierr); 208 ierr = KSPGetPC(H0ksp, &H0pc);CHKERRQ(ierr); 209 ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0pc, "tao_h0_");CHKERRQ(ierr); 210 211 ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr); 212 ierr = KSPSetUp(H0ksp);CHKERRQ(ierr); 213 } 214 215 PetscFunctionReturn(0); 216 } 217 218 /* ---------------------------------------------------------- */ 219 static PetscErrorCode TaoDestroy_LMVM(Tao tao) 220 { 221 TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 222 PetscErrorCode ierr; 223 PetscBool recycle; 224 225 PetscFunctionBegin; 226 if (lmP->M) { 227 ierr = MatLMVMGetRecycleFlag(lmP->M, &recycle); CHKERRQ(ierr); 228 if (recycle) { 229 ierr = PetscInfo(tao, "WARNING: TaoDestroy() called when LMVM recycling is enabled!\n"); CHKERRQ(ierr); 230 } 231 } 232 233 if (tao->setupcalled) { 234 ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 235 ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 236 ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 237 ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 238 } 239 240 if (lmP->H0) { 241 ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr); 242 } 243 244 ierr = PetscFree(tao->data);CHKERRQ(ierr); 245 246 PetscFunctionReturn(0); 247 } 248 249 /*------------------------------------------------------------*/ 250 static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao) 251 { 252 PetscErrorCode ierr; 253 254 PetscFunctionBegin; 255 ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr); 256 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 257 ierr = PetscOptionsTail();CHKERRQ(ierr); 258 PetscFunctionReturn(0); 259 } 260 261 /*------------------------------------------------------------*/ 262 static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) 263 { 264 TAO_LMVM *lm = (TAO_LMVM *)tao->data; 265 PetscBool isascii, recycle; 266 PetscInt recycled_its; 267 PetscErrorCode ierr; 268 269 PetscFunctionBegin; 270 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 271 if (isascii) { 272 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 273 ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);CHKERRQ(ierr); 274 ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);CHKERRQ(ierr); 275 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 276 ierr = MatLMVMGetRecycleFlag(lm->M, &recycle);CHKERRQ(ierr); 277 if (recycle) { 278 ierr = PetscViewerASCIIPrintf(viewer, "Recycle: on\n");CHKERRQ(ierr); 279 recycled_its = lm->bfgs + lm->sgrad + lm->grad; 280 ierr = PetscViewerASCIIPrintf(viewer, "Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr); 281 } 282 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 283 } 284 PetscFunctionReturn(0); 285 } 286 287 /* ---------------------------------------------------------- */ 288 289 /*MC 290 TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton 291 optimization solver for unconstrained minimization. It solves 292 the Newton step 293 Hkdk = - gk 294 295 using an approximation Bk in place of Hk, where Bk is composed using 296 the BFGS update formula. A More-Thuente line search is then used 297 to computed the steplength in the dk direction 298 Options Database Keys: 299 + -tao_lmm_vectors - number of vectors to use for approximation 300 . -tao_lmm_scale_type - "none","scalar","broyden" 301 . -tao_lmm_limit_type - "none","average","relative","absolute" 302 . -tao_lmm_rescale_type - "none","scalar","gl" 303 . -tao_lmm_limit_mu - mu limiting factor 304 . -tao_lmm_limit_nu - nu limiting factor 305 . -tao_lmm_delta_min - minimum delta value 306 . -tao_lmm_delta_max - maximum delta value 307 . -tao_lmm_broyden_phi - phi factor for Broyden scaling 308 . -tao_lmm_scalar_alpha - alpha factor for scalar scaling 309 . -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal 310 . -tao_lmm_rescale_beta - beta factor for rescaling diagonal 311 . -tao_lmm_scalar_history - amount of history for scalar scaling 312 . -tao_lmm_rescale_history - amount of history for rescaling diagonal 313 . -tao_lmm_eps - rejection tolerance 314 - -tao_lmm_recycle - enable recycling LMVM updates between TaoSolve() calls 315 316 Level: beginner 317 M*/ 318 319 PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao) 320 { 321 TAO_LMVM *lmP; 322 const char *morethuente_type = TAOLINESEARCHMT; 323 PetscErrorCode ierr; 324 325 PetscFunctionBegin; 326 tao->ops->setup = TaoSetUp_LMVM; 327 tao->ops->solve = TaoSolve_LMVM; 328 tao->ops->view = TaoView_LMVM; 329 tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 330 tao->ops->destroy = TaoDestroy_LMVM; 331 332 ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 333 lmP->D = 0; 334 lmP->M = 0; 335 lmP->Xold = 0; 336 lmP->Gold = 0; 337 lmP->H0 = NULL; 338 339 tao->data = (void*)lmP; 340 /* Override default settings (unless already changed) */ 341 if (!tao->max_it_changed) tao->max_it = 2000; 342 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 343 344 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 345 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 346 ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 347 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 348 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 349 PetscFunctionReturn(0); 350 } 351