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