1 #include <petsctaolinesearch.h> 2 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 3 4 #define LMVM_STEP_BFGS 0 5 #define LMVM_STEP_GRAD 1 6 7 static PetscErrorCode TaoSolve_LMVM(Tao tao) 8 { 9 TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 10 PetscReal f, fold, gdx, gnorm; 11 PetscReal step = 1.0; 12 PetscErrorCode ierr; 13 PetscInt stepType = LMVM_STEP_GRAD, nupdates; 14 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 15 16 PetscFunctionBegin; 17 18 if (tao->XL || tao->XU || tao->ops->computebounds) { 19 ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr); 20 } 21 22 /* Check convergence criteria */ 23 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 24 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 25 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,step);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 (!lmP->recycle) { 36 lmP->bfgs = 0; 37 lmP->grad = 0; 38 ierr = MatLMVMReset(lmP->M, PETSC_FALSE); CHKERRQ(ierr); 39 } 40 41 /* Have not converged; continue with Newton method */ 42 while (tao->reason == TAO_CONTINUE_ITERATING) { 43 /* Compute direction */ 44 if (lmP->H0) { 45 ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 46 stepType = LMVM_STEP_BFGS; 47 } 48 ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 49 ierr = MatSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 50 ierr = MatLMVMGetUpdateCount(lmP->M, &nupdates); CHKERRQ(ierr); 51 if (nupdates > 0) stepType = LMVM_STEP_BFGS; 52 53 /* Check for success (descent direction) */ 54 ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr); 55 if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 56 /* Step is not descent or direction produced not a number 57 We can assert bfgsUpdates > 1 in this case because 58 the first solve produces the scaled gradient direction, 59 which is guaranteed to be descent 60 61 Use steepest descent direction (scaled) 62 */ 63 64 ierr = MatLMVMResetJ0(lmP->M);CHKERRQ(ierr); 65 ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 66 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 67 ierr = MatSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr); 68 69 /* On a reset, the direction cannot be not a number; it is a 70 scaled gradient step. No need to check for this condition. */ 71 stepType = LMVM_STEP_GRAD; 72 } 73 ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 74 75 /* Perform the linesearch */ 76 fold = f; 77 ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 78 ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 79 80 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr); 81 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 82 83 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) { 84 /* Reset factors and use scaled gradient step */ 85 f = fold; 86 ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 87 ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 88 89 /* Failed to obtain acceptable iterate with BFGS step */ 90 /* Attempt to use the scaled gradient direction */ 91 92 ierr = MatLMVMResetJ0(lmP->M);CHKERRQ(ierr); 93 ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 94 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 95 ierr = MatSolve(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 96 97 /* On a reset, the direction cannot be not a number; it is a 98 scaled gradient step. No need to check for this condition. */ 99 stepType = LMVM_STEP_GRAD; 100 ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 101 102 /* Perform the linesearch */ 103 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr); 104 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 105 } 106 107 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 108 /* Failed to find an improving point */ 109 f = fold; 110 ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 111 ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 112 step = 0.0; 113 tao->reason = TAO_DIVERGED_LS_FAILURE; 114 } else { 115 /* LS found valid step, so tally up step type */ 116 switch (stepType) { 117 case LMVM_STEP_BFGS: 118 ++lmP->bfgs; 119 break; 120 case LMVM_STEP_GRAD: 121 ++lmP->grad; 122 break; 123 default: 124 break; 125 } 126 /* Compute new gradient norm */ 127 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 128 } 129 130 /* Check convergence */ 131 tao->niter++; 132 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 133 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 134 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 135 } 136 PetscFunctionReturn(0); 137 } 138 139 static PetscErrorCode TaoSetUp_LMVM(Tao tao) 140 { 141 TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 142 PetscInt n,N; 143 PetscErrorCode ierr; 144 145 PetscFunctionBegin; 146 /* Existence of tao->solution checked in TaoSetUp() */ 147 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 148 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 149 if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 150 if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 151 if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 152 153 /* Create matrix for the limited memory approximation */ 154 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 155 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 156 ierr = MatSetSizes(lmP->M, n, n, N, N);CHKERRQ(ierr); 157 ierr = MatLMVMAllocate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 158 159 /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 160 if (lmP->H0) { 161 ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 162 } 163 164 PetscFunctionReturn(0); 165 } 166 167 /* ---------------------------------------------------------- */ 168 static PetscErrorCode TaoDestroy_LMVM(Tao tao) 169 { 170 TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 171 PetscErrorCode ierr; 172 173 PetscFunctionBegin; 174 if (tao->setupcalled) { 175 ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 176 ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 177 ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 178 } 179 ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 180 if (lmP->H0) { 181 ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr); 182 } 183 ierr = PetscFree(tao->data);CHKERRQ(ierr); 184 185 PetscFunctionReturn(0); 186 } 187 188 /*------------------------------------------------------------*/ 189 static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao) 190 { 191 TAO_LMVM *lm = (TAO_LMVM *)tao->data; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr); 196 ierr = PetscOptionsBool("-tao_lmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",lm->recycle,&lm->recycle,NULL);CHKERRQ(ierr); 197 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 198 ierr = MatSetFromOptions(lm->M);CHKERRQ(ierr); 199 ierr = PetscOptionsTail();CHKERRQ(ierr); 200 PetscFunctionReturn(0); 201 } 202 203 /*------------------------------------------------------------*/ 204 static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) 205 { 206 TAO_LMVM *lm = (TAO_LMVM *)tao->data; 207 PetscBool isascii; 208 PetscInt recycled_its; 209 PetscErrorCode ierr; 210 211 PetscFunctionBegin; 212 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 213 if (isascii) { 214 ierr = PetscViewerASCIIPrintf(viewer, " Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 215 if (lm->recycle) { 216 ierr = PetscViewerASCIIPrintf(viewer, " Recycle: on\n");CHKERRQ(ierr); 217 recycled_its = lm->bfgs + lm->grad; 218 ierr = PetscViewerASCIIPrintf(viewer, " Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr); 219 } 220 } 221 PetscFunctionReturn(0); 222 } 223 224 /* ---------------------------------------------------------- */ 225 226 /*MC 227 TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton 228 optimization solver for unconstrained minimization. It solves 229 the Newton step 230 Hkdk = - gk 231 232 using an approximation Bk in place of Hk, where Bk is composed using 233 the BFGS update formula. A More-Thuente line search is then used 234 to computed the steplength in the dk direction 235 Options Database Keys: 236 . -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls 237 238 Level: beginner 239 M*/ 240 241 PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao) 242 { 243 TAO_LMVM *lmP; 244 const char *morethuente_type = TAOLINESEARCHMT; 245 PetscErrorCode ierr; 246 247 PetscFunctionBegin; 248 tao->ops->setup = TaoSetUp_LMVM; 249 tao->ops->solve = TaoSolve_LMVM; 250 tao->ops->view = TaoView_LMVM; 251 tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 252 tao->ops->destroy = TaoDestroy_LMVM; 253 254 ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 255 lmP->D = 0; 256 lmP->M = 0; 257 lmP->Xold = 0; 258 lmP->Gold = 0; 259 lmP->H0 = NULL; 260 lmP->recycle = PETSC_FALSE; 261 262 tao->data = (void*)lmP; 263 /* Override default settings (unless already changed) */ 264 if (!tao->max_it_changed) tao->max_it = 2000; 265 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 266 267 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 268 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 269 ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 270 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 271 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 272 273 ierr = KSPInitializePackage();CHKERRQ(ierr); 274 ierr = MatCreate(((PetscObject)tao)->comm, &lmP->M);CHKERRQ(ierr); 275 ierr = PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1);CHKERRQ(ierr); 276 ierr = MatSetType(lmP->M, MATLMVMBFGS);CHKERRQ(ierr); 277 ierr = MatSetOptionsPrefix(lmP->M, "tao_lmvm_");CHKERRQ(ierr); 278 PetscFunctionReturn(0); 279 } 280