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