1 #include <petsctaolinesearch.h> 2 #include <../src/tao/matrix/lmvmmat.h> 3 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 4 #include <../src/tao/bound/impls/blmvm/blmvm.h> 5 6 /*------------------------------------------------------------*/ 7 static PetscErrorCode TaoSolve_BLMVM(Tao tao) 8 { 9 PetscErrorCode ierr; 10 TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 11 TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 12 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 13 PetscReal f, fold, gdx, gnorm; 14 PetscReal stepsize = 1.0,delta; 15 16 PetscFunctionBegin; 17 /* Project initial point onto bounds */ 18 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 19 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 20 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 21 22 23 /* Check convergence criteria */ 24 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr); 25 ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 26 27 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 28 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 29 30 ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr); 31 if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 32 33 /* Set initial scaling for the function */ 34 if (f != 0.0) { 35 delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm); 36 } else { 37 delta = 2.0 / (gnorm*gnorm); 38 } 39 ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr); 40 ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr); 41 42 /* Set counter for gradient/reset steps */ 43 blmP->grad = 0; 44 blmP->reset = 0; 45 46 /* Have not converged; continue with Newton method */ 47 while (reason == TAO_CONTINUE_ITERATING) { 48 /* Compute direction */ 49 ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 50 ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 51 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 52 53 /* Check for success (descent direction) */ 54 ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr); 55 if (gdx <= 0) { 56 /* Step is not descent or solve was not successful 57 Use steepest descent direction (scaled) */ 58 ++blmP->grad; 59 60 if (f != 0.0) { 61 delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm); 62 } else { 63 delta = 2.0 / (gnorm*gnorm); 64 } 65 ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr); 66 ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr); 67 ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 68 ierr = MatLMVMSolve(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 if (f != 0.0) { 90 delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm); 91 } else { 92 delta = 2.0/ (gnorm*gnorm); 93 } 94 ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr); 95 ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr); 96 ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 97 ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 98 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 99 100 /* This may be incorrect; linesearch has values for stepmax and stepmin 101 that should be reset. */ 102 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 103 ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr); 104 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 105 106 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 107 tao->reason = TAO_DIVERGED_LS_FAILURE; 108 break; 109 } 110 } 111 112 /* Check for converged */ 113 ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 114 ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 115 116 117 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number"); 118 tao->niter++; 119 ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr); 120 } 121 PetscFunctionReturn(0); 122 } 123 124 static PetscErrorCode TaoSetup_BLMVM(Tao tao) 125 { 126 TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 127 PetscInt n,N; 128 PetscErrorCode ierr; 129 KSP H0ksp; 130 131 PetscFunctionBegin; 132 /* Existence of tao->solution checked in TaoSetup() */ 133 ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr); 134 ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr); 135 ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr); 136 137 if (!tao->stepdirection) { 138 ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 139 } 140 if (!tao->gradient) { 141 ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 142 } 143 if (!tao->XL) { 144 ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr); 145 ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr); 146 } 147 if (!tao->XU) { 148 ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr); 149 ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr); 150 } 151 /* Create matrix for the limited memory approximation */ 152 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 153 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 154 ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);CHKERRQ(ierr); 155 ierr = MatLMVMAllocateVectors(blmP->M,tao->solution);CHKERRQ(ierr); 156 157 /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 158 if (blmP->H0) { 159 const char *prefix; 160 PC H0pc; 161 162 ierr = MatLMVMSetH0(blmP->M, blmP->H0);CHKERRQ(ierr); 163 ierr = MatLMVMGetH0KSP(blmP->M, &H0ksp);CHKERRQ(ierr); 164 165 ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr); 166 ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr); 167 ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0ksp, "tao_h0_");CHKERRQ(ierr); 168 ierr = KSPGetPC(H0ksp, &H0pc);CHKERRQ(ierr); 169 ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0pc, "tao_h0_");CHKERRQ(ierr); 170 171 ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr); 172 ierr = KSPSetUp(H0ksp);CHKERRQ(ierr); 173 } 174 175 PetscFunctionReturn(0); 176 } 177 178 /* ---------------------------------------------------------- */ 179 static PetscErrorCode TaoDestroy_BLMVM(Tao tao) 180 { 181 TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 182 PetscErrorCode ierr; 183 184 PetscFunctionBegin; 185 if (tao->setupcalled) { 186 ierr = MatDestroy(&blmP->M);CHKERRQ(ierr); 187 ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr); 188 ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr); 189 ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr); 190 } 191 192 if (blmP->H0) { 193 PetscObjectDereference((PetscObject)blmP->H0); 194 } 195 196 ierr = PetscFree(tao->data);CHKERRQ(ierr); 197 PetscFunctionReturn(0); 198 } 199 200 /*------------------------------------------------------------*/ 201 static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao) 202 { 203 PetscErrorCode ierr; 204 205 PetscFunctionBegin; 206 ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr); 207 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 208 ierr = PetscOptionsTail();CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 213 /*------------------------------------------------------------*/ 214 static int TaoView_BLMVM(Tao tao, PetscViewer viewer) 215 { 216 TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 217 PetscBool isascii; 218 PetscErrorCode ierr; 219 220 PetscFunctionBegin; 221 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 222 if (isascii) { 223 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 224 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr); 225 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 226 } 227 PetscFunctionReturn(0); 228 } 229 230 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 231 { 232 TAO_BLMVM *blm = (TAO_BLMVM *) tao->data; 233 PetscErrorCode ierr; 234 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 237 PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 238 PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 239 if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 240 241 ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr); 242 ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr); 243 ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 244 ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 245 246 ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr); 247 ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr); 248 ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr); 249 PetscFunctionReturn(0); 250 } 251 252 /* ---------------------------------------------------------- */ 253 /*MC 254 TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 255 for nonlinear minimization with bound constraints. It is an extension 256 of TAOLMVM 257 258 Options Database Keys: 259 + -tao_lmm_vectors - number of vectors to use for approximation 260 . -tao_lmm_scale_type - "none","scalar","broyden" 261 . -tao_lmm_limit_type - "none","average","relative","absolute" 262 . -tao_lmm_rescale_type - "none","scalar","gl" 263 . -tao_lmm_limit_mu - mu limiting factor 264 . -tao_lmm_limit_nu - nu limiting factor 265 . -tao_lmm_delta_min - minimum delta value 266 . -tao_lmm_delta_max - maximum delta value 267 . -tao_lmm_broyden_phi - phi factor for Broyden scaling 268 . -tao_lmm_scalar_alpha - alpha factor for scalar scaling 269 . -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal 270 . -tao_lmm_rescale_beta - beta factor for rescaling diagonal 271 . -tao_lmm_scalar_history - amount of history for scalar scaling 272 . -tao_lmm_rescale_history - amount of history for rescaling diagonal 273 - -tao_lmm_eps - rejection tolerance 274 275 Level: beginner 276 M*/ 277 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 278 { 279 TAO_BLMVM *blmP; 280 const char *morethuente_type = TAOLINESEARCHMT; 281 PetscErrorCode ierr; 282 283 PetscFunctionBegin; 284 tao->ops->setup = TaoSetup_BLMVM; 285 tao->ops->solve = TaoSolve_BLMVM; 286 tao->ops->view = TaoView_BLMVM; 287 tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 288 tao->ops->destroy = TaoDestroy_BLMVM; 289 tao->ops->computedual = TaoComputeDual_BLMVM; 290 291 ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr); 292 blmP->H0 = NULL; 293 tao->data = (void*)blmP; 294 295 /* Override default settings (unless already changed) */ 296 if (!tao->max_it_changed) tao->max_it = 2000; 297 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 298 299 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 300 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 301 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 302 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) 307 { 308 TAO_LMVM *lmP; 309 TAO_BLMVM *blmP; 310 const TaoType type; 311 PetscBool is_lmvm, is_blmvm; 312 313 PetscErrorCode ierr; 314 315 ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 316 ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 317 ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 318 319 if (is_lmvm) { 320 lmP = (TAO_LMVM *)tao->data; 321 ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 322 lmP->H0 = H0; 323 } else if (is_blmvm) { 324 blmP = (TAO_BLMVM *)tao->data; 325 ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 326 blmP->H0 = H0; 327 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 328 329 PetscFunctionReturn(0); 330 } 331 332 PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) 333 { 334 TAO_LMVM *lmP; 335 TAO_BLMVM *blmP; 336 const TaoType type; 337 PetscBool is_lmvm, is_blmvm; 338 Mat M; 339 340 PetscErrorCode ierr; 341 342 ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 343 ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 344 ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 345 346 if (is_lmvm) { 347 lmP = (TAO_LMVM *)tao->data; 348 M = lmP->M; 349 } else if (is_blmvm) { 350 blmP = (TAO_BLMVM *)tao->data; 351 M = blmP->M; 352 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 353 354 ierr = MatLMVMGetH0(M, H0);CHKERRQ(ierr); 355 PetscFunctionReturn(0); 356 } 357 358 PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) 359 { 360 TAO_LMVM *lmP; 361 TAO_BLMVM *blmP; 362 const TaoType type; 363 PetscBool is_lmvm, is_blmvm; 364 Mat M; 365 PetscErrorCode ierr; 366 367 ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 368 ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 369 ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 370 371 if (is_lmvm) { 372 lmP = (TAO_LMVM *)tao->data; 373 M = lmP->M; 374 } else if (is_blmvm) { 375 blmP = (TAO_BLMVM *)tao->data; 376 M = blmP->M; 377 } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 378 379 ierr = MatLMVMGetH0KSP(M, ksp);CHKERRQ(ierr); 380 PetscFunctionReturn(0); 381 } 382