1 #include <petscksp.h> 2 #include <../src/tao/quadratic/impls/gpcg/gpcg.h> /*I "gpcg.h" I*/ 3 4 static PetscErrorCode GPCGGradProjections(Tao tao); 5 static PetscErrorCode GPCGObjectiveAndGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *); 6 7 /*------------------------------------------------------------*/ 8 static PetscErrorCode TaoDestroy_GPCG(Tao tao) 9 { 10 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 11 12 /* Free allocated memory in GPCG structure */ 13 PetscFunctionBegin; 14 PetscCall(VecDestroy(&gpcg->B)); 15 PetscCall(VecDestroy(&gpcg->Work)); 16 PetscCall(VecDestroy(&gpcg->X_New)); 17 PetscCall(VecDestroy(&gpcg->G_New)); 18 PetscCall(VecDestroy(&gpcg->DXFree)); 19 PetscCall(VecDestroy(&gpcg->R)); 20 PetscCall(VecDestroy(&gpcg->PG)); 21 PetscCall(MatDestroy(&gpcg->Hsub)); 22 PetscCall(MatDestroy(&gpcg->Hsub_pre)); 23 PetscCall(ISDestroy(&gpcg->Free_Local)); 24 PetscCall(KSPDestroy(&tao->ksp)); 25 PetscCall(PetscFree(tao->data)); 26 PetscFunctionReturn(PETSC_SUCCESS); 27 } 28 29 /*------------------------------------------------------------*/ 30 static PetscErrorCode TaoSetFromOptions_GPCG(Tao tao, PetscOptionItems PetscOptionsObject) 31 { 32 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 33 PetscBool flg; 34 35 PetscFunctionBegin; 36 PetscOptionsHeadBegin(PetscOptionsObject, "Gradient Projection, Conjugate Gradient method for bound constrained optimization"); 37 PetscCall(PetscOptionsInt("-tao_gpcg_maxpgits", "maximum number of gradient projections per GPCG iterate", NULL, gpcg->maxgpits, &gpcg->maxgpits, &flg)); 38 PetscOptionsHeadEnd(); 39 PetscCall(KSPSetFromOptions(tao->ksp)); 40 PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 41 PetscFunctionReturn(PETSC_SUCCESS); 42 } 43 44 /*------------------------------------------------------------*/ 45 static PetscErrorCode TaoView_GPCG(Tao tao, PetscViewer viewer) 46 { 47 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 48 PetscBool isascii; 49 50 PetscFunctionBegin; 51 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 52 if (isascii) { 53 PetscCall(PetscViewerASCIIPrintf(viewer, "Total PG its: %" PetscInt_FMT ",", gpcg->total_gp_its)); 54 PetscCall(PetscViewerASCIIPrintf(viewer, "PG tolerance: %g \n", (double)gpcg->pg_ftol)); 55 } 56 PetscCall(TaoLineSearchView(tao->linesearch, viewer)); 57 PetscFunctionReturn(PETSC_SUCCESS); 58 } 59 60 /* GPCGObjectiveAndGradient() 61 Compute f=0.5 * x'Hx + b'x + c 62 g=Hx + b 63 */ 64 static PetscErrorCode GPCGObjectiveAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *tptr) 65 { 66 Tao tao = (Tao)tptr; 67 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 68 PetscReal f1, f2; 69 70 PetscFunctionBegin; 71 PetscCall(MatMult(tao->hessian, X, G)); 72 PetscCall(VecDot(G, X, &f1)); 73 PetscCall(VecDot(gpcg->B, X, &f2)); 74 PetscCall(VecAXPY(G, 1.0, gpcg->B)); 75 *f = f1 / 2.0 + f2 + gpcg->c; 76 PetscFunctionReturn(PETSC_SUCCESS); 77 } 78 79 /* ---------------------------------------------------------- */ 80 static PetscErrorCode TaoSetup_GPCG(Tao tao) 81 { 82 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 83 84 PetscFunctionBegin; 85 /* Allocate some arrays */ 86 if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 87 if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 88 89 PetscCall(VecDuplicate(tao->solution, &gpcg->B)); 90 PetscCall(VecDuplicate(tao->solution, &gpcg->Work)); 91 PetscCall(VecDuplicate(tao->solution, &gpcg->X_New)); 92 PetscCall(VecDuplicate(tao->solution, &gpcg->G_New)); 93 PetscCall(VecDuplicate(tao->solution, &gpcg->DXFree)); 94 PetscCall(VecDuplicate(tao->solution, &gpcg->R)); 95 PetscCall(VecDuplicate(tao->solution, &gpcg->PG)); 96 /* 97 if (gpcg->ksp_type == GPCG_KSP_NASH) { 98 PetscCall(KSPSetType(tao->ksp,KSPNASH)); 99 } else if (gpcg->ksp_type == GPCG_KSP_STCG) { 100 PetscCall(KSPSetType(tao->ksp,KSPSTCG)); 101 } else { 102 PetscCall(KSPSetType(tao->ksp,KSPGLTR)); 103 } 104 if (tao->ksp->ops->setfromoptions) (*tao->ksp->ops->setfromoptions)(tao->ksp); 105 106 } 107 */ 108 PetscFunctionReturn(PETSC_SUCCESS); 109 } 110 111 static PetscErrorCode TaoSolve_GPCG(Tao tao) 112 { 113 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 114 PetscInt its; 115 PetscReal actred, f, f_new, gnorm, gdx, stepsize, xtb; 116 PetscReal xtHx; 117 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 118 119 PetscFunctionBegin; 120 PetscCall(TaoComputeVariableBounds(tao)); 121 PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 122 PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU)); 123 124 /* Using f = .5*x'Hx + x'b + c and g=Hx + b, compute b,c */ 125 PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre)); 126 PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 127 PetscCall(VecCopy(tao->gradient, gpcg->B)); 128 PetscCall(MatMult(tao->hessian, tao->solution, gpcg->Work)); 129 PetscCall(VecDot(gpcg->Work, tao->solution, &xtHx)); 130 PetscCall(VecAXPY(gpcg->B, -1.0, gpcg->Work)); 131 PetscCall(VecDot(gpcg->B, tao->solution, &xtb)); 132 gpcg->c = f - xtHx / 2.0 - xtb; 133 if (gpcg->Free_Local) PetscCall(ISDestroy(&gpcg->Free_Local)); 134 PetscCall(VecWhichInactive(tao->XL, tao->solution, tao->gradient, tao->XU, PETSC_TRUE, &gpcg->Free_Local)); 135 136 /* Project the gradient and calculate the norm */ 137 PetscCall(VecCopy(tao->gradient, gpcg->G_New)); 138 PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->PG)); 139 PetscCall(VecNorm(gpcg->PG, NORM_2, &gpcg->gnorm)); 140 tao->step = 1.0; 141 gpcg->f = f; 142 143 /* Check Stopping Condition */ 144 tao->reason = TAO_CONTINUE_ITERATING; 145 PetscCall(TaoLogConvergenceHistory(tao, f, gpcg->gnorm, 0.0, tao->ksp_its)); 146 PetscCall(TaoMonitor(tao, tao->niter, f, gpcg->gnorm, 0.0, tao->step)); 147 PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 148 149 while (tao->reason == TAO_CONTINUE_ITERATING) { 150 /* Call general purpose update function */ 151 PetscTryTypeMethod(tao, update, tao->niter, tao->user_update); 152 tao->ksp_its = 0; 153 154 PetscCall(GPCGGradProjections(tao)); 155 PetscCall(ISGetSize(gpcg->Free_Local, &gpcg->n_free)); 156 157 f = gpcg->f; 158 gnorm = gpcg->gnorm; 159 160 PetscCall(KSPReset(tao->ksp)); 161 162 if (gpcg->n_free > 0) { 163 /* Create a reduced linear system */ 164 PetscCall(VecDestroy(&gpcg->R)); 165 PetscCall(VecDestroy(&gpcg->DXFree)); 166 PetscCall(TaoVecGetSubVec(tao->gradient, gpcg->Free_Local, tao->subset_type, 0.0, &gpcg->R)); 167 PetscCall(VecScale(gpcg->R, -1.0)); 168 PetscCall(TaoVecGetSubVec(tao->stepdirection, gpcg->Free_Local, tao->subset_type, 0.0, &gpcg->DXFree)); 169 PetscCall(VecSet(gpcg->DXFree, 0.0)); 170 171 PetscCall(TaoMatGetSubMat(tao->hessian, gpcg->Free_Local, gpcg->Work, tao->subset_type, &gpcg->Hsub)); 172 173 if (tao->hessian_pre == tao->hessian) { 174 PetscCall(MatDestroy(&gpcg->Hsub_pre)); 175 PetscCall(PetscObjectReference((PetscObject)gpcg->Hsub)); 176 gpcg->Hsub_pre = gpcg->Hsub; 177 } else { 178 PetscCall(TaoMatGetSubMat(tao->hessian, gpcg->Free_Local, gpcg->Work, tao->subset_type, &gpcg->Hsub_pre)); 179 } 180 181 PetscCall(KSPReset(tao->ksp)); 182 PetscCall(KSPSetOperators(tao->ksp, gpcg->Hsub, gpcg->Hsub_pre)); 183 184 PetscCall(KSPSolve(tao->ksp, gpcg->R, gpcg->DXFree)); 185 PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 186 tao->ksp_its += its; 187 tao->ksp_tot_its += its; 188 PetscCall(VecSet(tao->stepdirection, 0.0)); 189 PetscCall(VecISAXPY(tao->stepdirection, gpcg->Free_Local, 1.0, gpcg->DXFree)); 190 191 PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx)); 192 PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 193 f_new = f; 194 PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, &stepsize, &ls_status)); 195 196 actred = f_new - f; 197 198 /* Evaluate the function and gradient at the new point */ 199 PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->PG)); 200 PetscCall(VecNorm(gpcg->PG, NORM_2, &gnorm)); 201 f = f_new; 202 PetscCall(ISDestroy(&gpcg->Free_Local)); 203 PetscCall(VecWhichInactive(tao->XL, tao->solution, tao->gradient, tao->XU, PETSC_TRUE, &gpcg->Free_Local)); 204 } else { 205 actred = 0; 206 gpcg->step = 1.0; 207 /* if there were no free variables, no cg method */ 208 } 209 210 tao->niter++; 211 gpcg->f = f; 212 gpcg->gnorm = gnorm; 213 gpcg->actred = actred; 214 PetscCall(TaoLogConvergenceHistory(tao, f, gpcg->gnorm, 0.0, tao->ksp_its)); 215 PetscCall(TaoMonitor(tao, tao->niter, f, gpcg->gnorm, 0.0, tao->step)); 216 PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 217 if (tao->reason != TAO_CONTINUE_ITERATING) break; 218 } /* END MAIN LOOP */ 219 PetscFunctionReturn(PETSC_SUCCESS); 220 } 221 222 static PetscErrorCode GPCGGradProjections(Tao tao) 223 { 224 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 225 PetscInt i; 226 PetscReal actred = -1.0, actred_max = 0.0, gAg, gtg = gpcg->gnorm, alpha; 227 PetscReal f_new, gdx, stepsize; 228 Vec DX = tao->stepdirection, XL = tao->XL, XU = tao->XU, Work = gpcg->Work; 229 Vec X = tao->solution, G = tao->gradient; 230 TaoLineSearchConvergedReason lsflag = TAOLINESEARCH_CONTINUE_ITERATING; 231 232 /* 233 The free, active, and binding variables should be already identified 234 */ 235 PetscFunctionBegin; 236 for (i = 0; i < gpcg->maxgpits; i++) { 237 if (-actred <= (gpcg->pg_ftol) * actred_max) break; 238 PetscCall(VecBoundGradientProjection(G, X, XL, XU, DX)); 239 PetscCall(VecScale(DX, -1.0)); 240 PetscCall(VecDot(DX, G, &gdx)); 241 242 PetscCall(MatMult(tao->hessian, DX, Work)); 243 PetscCall(VecDot(DX, Work, &gAg)); 244 245 gpcg->gp_iterates++; 246 gpcg->total_gp_its++; 247 248 gtg = -gdx; 249 if (PetscAbsReal(gAg) == 0.0) { 250 alpha = 1.0; 251 } else { 252 alpha = PetscAbsReal(gtg / gAg); 253 } 254 PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, alpha)); 255 f_new = gpcg->f; 256 PetscCall(TaoLineSearchApply(tao->linesearch, X, &f_new, G, DX, &stepsize, &lsflag)); 257 258 /* Update the iterate */ 259 actred = f_new - gpcg->f; 260 actred_max = PetscMax(actred_max, -(f_new - gpcg->f)); 261 gpcg->f = f_new; 262 PetscCall(ISDestroy(&gpcg->Free_Local)); 263 PetscCall(VecWhichInactive(XL, X, tao->gradient, XU, PETSC_TRUE, &gpcg->Free_Local)); 264 } 265 266 gpcg->gnorm = gtg; 267 PetscFunctionReturn(PETSC_SUCCESS); 268 } /* End gradient projections */ 269 270 static PetscErrorCode TaoComputeDual_GPCG(Tao tao, Vec DXL, Vec DXU) 271 { 272 TAO_GPCG *gpcg = (TAO_GPCG *)tao->data; 273 274 PetscFunctionBegin; 275 PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->Work)); 276 PetscCall(VecCopy(gpcg->Work, DXL)); 277 PetscCall(VecAXPY(DXL, -1.0, tao->gradient)); 278 PetscCall(VecSet(DXU, 0.0)); 279 PetscCall(VecPointwiseMax(DXL, DXL, DXU)); 280 281 PetscCall(VecCopy(tao->gradient, DXU)); 282 PetscCall(VecAXPY(DXU, -1.0, gpcg->Work)); 283 PetscCall(VecSet(gpcg->Work, 0.0)); 284 PetscCall(VecPointwiseMin(DXU, gpcg->Work, DXU)); 285 PetscFunctionReturn(PETSC_SUCCESS); 286 } 287 288 /*------------------------------------------------------------*/ 289 /*MC 290 TAOGPCG - gradient projected conjugate gradient algorithm is an active-set 291 conjugate-gradient based method for bound-constrained minimization 292 293 Options Database Keys: 294 + -tao_gpcg_maxpgits - maximum number of gradient projections for GPCG iterate 295 - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 296 297 Level: beginner 298 M*/ 299 PETSC_EXTERN PetscErrorCode TaoCreate_GPCG(Tao tao) 300 { 301 TAO_GPCG *gpcg; 302 303 PetscFunctionBegin; 304 tao->ops->setup = TaoSetup_GPCG; 305 tao->ops->solve = TaoSolve_GPCG; 306 tao->ops->view = TaoView_GPCG; 307 tao->ops->setfromoptions = TaoSetFromOptions_GPCG; 308 tao->ops->destroy = TaoDestroy_GPCG; 309 tao->ops->computedual = TaoComputeDual_GPCG; 310 311 PetscCall(PetscNew(&gpcg)); 312 tao->data = (void *)gpcg; 313 314 /* Override default settings (unless already changed) */ 315 PetscCall(TaoParametersInitialize(tao)); 316 PetscObjectParameterSetDefault(tao, max_it, 500); 317 PetscObjectParameterSetDefault(tao, max_funcs, 100000); 318 PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1e-6 : 1e-12); 319 PetscObjectParameterSetDefault(tao, grtol, PetscDefined(USE_REAL_SINGLE) ? 1e-6 : 1e-12); 320 321 /* Initialize pointers and variables */ 322 gpcg->n = 0; 323 gpcg->maxgpits = 8; 324 gpcg->pg_ftol = 0.1; 325 326 gpcg->gp_iterates = 0; /* Cumulative number */ 327 gpcg->total_gp_its = 0; 328 329 /* Initialize pointers and variables */ 330 gpcg->n_bind = 0; 331 gpcg->n_free = 0; 332 gpcg->n_upper = 0; 333 gpcg->n_lower = 0; 334 gpcg->subset_type = TAO_SUBSET_MASK; 335 gpcg->Hsub = NULL; 336 gpcg->Hsub_pre = NULL; 337 338 PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 339 PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 340 PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 341 PetscCall(KSPSetType(tao->ksp, KSPNASH)); 342 343 PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 344 PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 345 PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHGPCG)); 346 PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, GPCGObjectiveAndGradient, tao)); 347 PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix)); 348 PetscFunctionReturn(PETSC_SUCCESS); 349 } 350