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