xref: /petsc/src/tao/quadratic/impls/gpcg/gpcg.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
174962610SAlp Dener #include <petscksp.h>
274962610SAlp Dener #include <../src/tao/quadratic/impls/gpcg/gpcg.h> /*I "gpcg.h" I*/
374962610SAlp Dener 
474962610SAlp Dener static PetscErrorCode GPCGGradProjections(Tao tao);
574962610SAlp Dener static PetscErrorCode GPCGObjectiveAndGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *);
674962610SAlp Dener 
TaoDestroy_GPCG(Tao tao)7d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_GPCG(Tao tao)
8d71ae5a4SJacob Faibussowitsch {
974962610SAlp Dener   TAO_GPCG *gpcg = (TAO_GPCG *)tao->data;
1074962610SAlp Dener 
1174962610SAlp Dener   /* Free allocated memory in GPCG structure */
1274962610SAlp Dener   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gpcg->B));
149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gpcg->Work));
159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gpcg->X_New));
169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gpcg->G_New));
179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gpcg->DXFree));
189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gpcg->R));
199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gpcg->PG));
209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&gpcg->Hsub));
219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&gpcg->Hsub_pre));
229566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&gpcg->Free_Local));
23a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
249566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2674962610SAlp Dener }
2774962610SAlp Dener 
TaoSetFromOptions_GPCG(Tao tao,PetscOptionItems PetscOptionsObject)28ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_GPCG(Tao tao, PetscOptionItems PetscOptionsObject)
29d71ae5a4SJacob Faibussowitsch {
3074962610SAlp Dener   TAO_GPCG *gpcg = (TAO_GPCG *)tao->data;
3174962610SAlp Dener   PetscBool flg;
3274962610SAlp Dener 
3374962610SAlp Dener   PetscFunctionBegin;
34d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Gradient Projection, Conjugate Gradient method for bound constrained optimization");
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_gpcg_maxpgits", "maximum number of gradient projections per GPCG iterate", NULL, gpcg->maxgpits, &gpcg->maxgpits, &flg));
36d0609cedSBarry Smith   PetscOptionsHeadEnd();
379566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
389566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4074962610SAlp Dener }
4174962610SAlp Dener 
TaoView_GPCG(Tao tao,PetscViewer viewer)42d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_GPCG(Tao tao, PetscViewer viewer)
43d71ae5a4SJacob Faibussowitsch {
4474962610SAlp Dener   TAO_GPCG *gpcg = (TAO_GPCG *)tao->data;
4574962610SAlp Dener   PetscBool isascii;
4674962610SAlp Dener 
4774962610SAlp Dener   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
4974962610SAlp Dener   if (isascii) {
5063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Total PG its: %" PetscInt_FMT ",", gpcg->total_gp_its));
519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "PG tolerance: %g \n", (double)gpcg->pg_ftol));
5274962610SAlp Dener   }
539566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchView(tao->linesearch, viewer));
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5574962610SAlp Dener }
5674962610SAlp Dener 
5774962610SAlp Dener /* GPCGObjectiveAndGradient()
5874962610SAlp Dener    Compute f=0.5 * x'Hx + b'x + c
5974962610SAlp Dener            g=Hx + b
6074962610SAlp Dener */
GPCGObjectiveAndGradient(TaoLineSearch ls,Vec X,PetscReal * f,Vec G,void * tptr)61d71ae5a4SJacob Faibussowitsch static PetscErrorCode GPCGObjectiveAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *tptr)
62d71ae5a4SJacob Faibussowitsch {
6374962610SAlp Dener   Tao       tao  = (Tao)tptr;
6474962610SAlp Dener   TAO_GPCG *gpcg = (TAO_GPCG *)tao->data;
6574962610SAlp Dener   PetscReal f1, f2;
6674962610SAlp Dener 
6774962610SAlp Dener   PetscFunctionBegin;
689566063dSJacob Faibussowitsch   PetscCall(MatMult(tao->hessian, X, G));
699566063dSJacob Faibussowitsch   PetscCall(VecDot(G, X, &f1));
709566063dSJacob Faibussowitsch   PetscCall(VecDot(gpcg->B, X, &f2));
719566063dSJacob Faibussowitsch   PetscCall(VecAXPY(G, 1.0, gpcg->B));
7274962610SAlp Dener   *f = f1 / 2.0 + f2 + gpcg->c;
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7474962610SAlp Dener }
7574962610SAlp Dener 
TaoSetup_GPCG(Tao tao)76d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_GPCG(Tao tao)
77d71ae5a4SJacob Faibussowitsch {
7874962610SAlp Dener   TAO_GPCG *gpcg = (TAO_GPCG *)tao->data;
7974962610SAlp Dener 
8074962610SAlp Dener   PetscFunctionBegin;
8174962610SAlp Dener   /* Allocate some arrays */
8248a46eb9SPierre Jolivet   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
8348a46eb9SPierre Jolivet   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
8474962610SAlp Dener 
859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &gpcg->B));
869566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &gpcg->Work));
879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &gpcg->X_New));
889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &gpcg->G_New));
899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &gpcg->DXFree));
909566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &gpcg->R));
919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &gpcg->PG));
9274962610SAlp Dener   /*
9374962610SAlp Dener     if (gpcg->ksp_type == GPCG_KSP_NASH) {
949566063dSJacob Faibussowitsch         PetscCall(KSPSetType(tao->ksp,KSPNASH));
9574962610SAlp Dener       } else if (gpcg->ksp_type == GPCG_KSP_STCG) {
969566063dSJacob Faibussowitsch         PetscCall(KSPSetType(tao->ksp,KSPSTCG));
9774962610SAlp Dener       } else {
989566063dSJacob Faibussowitsch         PetscCall(KSPSetType(tao->ksp,KSPGLTR));
9974962610SAlp Dener       }
100*3a7d0413SPierre Jolivet       if (tao->ksp->ops->setfromoptions) (*tao->ksp->ops->setfromoptions)(tao->ksp);
10174962610SAlp Dener 
10274962610SAlp Dener     }
10374962610SAlp Dener   */
1043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10574962610SAlp Dener }
10674962610SAlp Dener 
TaoSolve_GPCG(Tao tao)107d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_GPCG(Tao tao)
108d71ae5a4SJacob Faibussowitsch {
10974962610SAlp Dener   TAO_GPCG                    *gpcg = (TAO_GPCG *)tao->data;
11074962610SAlp Dener   PetscInt                     its;
11174962610SAlp Dener   PetscReal                    actred, f, f_new, gnorm, gdx, stepsize, xtb;
11274962610SAlp Dener   PetscReal                    xtHx;
11374962610SAlp Dener   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
11474962610SAlp Dener 
11574962610SAlp Dener   PetscFunctionBegin;
1169566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
1179566063dSJacob Faibussowitsch   PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution));
1189566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
11974962610SAlp Dener 
12074962610SAlp Dener   /* Using f = .5*x'Hx + x'b + c and g=Hx + b,  compute b,c */
1219566063dSJacob Faibussowitsch   PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
1229566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
1239566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, gpcg->B));
1249566063dSJacob Faibussowitsch   PetscCall(MatMult(tao->hessian, tao->solution, gpcg->Work));
1259566063dSJacob Faibussowitsch   PetscCall(VecDot(gpcg->Work, tao->solution, &xtHx));
1269566063dSJacob Faibussowitsch   PetscCall(VecAXPY(gpcg->B, -1.0, gpcg->Work));
1279566063dSJacob Faibussowitsch   PetscCall(VecDot(gpcg->B, tao->solution, &xtb));
12874962610SAlp Dener   gpcg->c = f - xtHx / 2.0 - xtb;
12948a46eb9SPierre Jolivet   if (gpcg->Free_Local) PetscCall(ISDestroy(&gpcg->Free_Local));
1309566063dSJacob Faibussowitsch   PetscCall(VecWhichInactive(tao->XL, tao->solution, tao->gradient, tao->XU, PETSC_TRUE, &gpcg->Free_Local));
13174962610SAlp Dener 
13274962610SAlp Dener   /* Project the gradient and calculate the norm */
1339566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, gpcg->G_New));
1349566063dSJacob Faibussowitsch   PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->PG));
1359566063dSJacob Faibussowitsch   PetscCall(VecNorm(gpcg->PG, NORM_2, &gpcg->gnorm));
13674962610SAlp Dener   tao->step = 1.0;
13774962610SAlp Dener   gpcg->f   = f;
13874962610SAlp Dener 
13974962610SAlp Dener   /* Check Stopping Condition      */
14074962610SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1419566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, gpcg->gnorm, 0.0, tao->ksp_its));
1429566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, f, gpcg->gnorm, 0.0, tao->step));
143dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
14474962610SAlp Dener 
14574962610SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
146e1e80dc8SAlp Dener     /* Call general purpose update function */
147dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
14874962610SAlp Dener     tao->ksp_its = 0;
14974962610SAlp Dener 
1509566063dSJacob Faibussowitsch     PetscCall(GPCGGradProjections(tao));
1519566063dSJacob Faibussowitsch     PetscCall(ISGetSize(gpcg->Free_Local, &gpcg->n_free));
15274962610SAlp Dener 
1539371c9d4SSatish Balay     f     = gpcg->f;
1549371c9d4SSatish Balay     gnorm = gpcg->gnorm;
15574962610SAlp Dener 
1569566063dSJacob Faibussowitsch     PetscCall(KSPReset(tao->ksp));
15774962610SAlp Dener 
15874962610SAlp Dener     if (gpcg->n_free > 0) {
15974962610SAlp Dener       /* Create a reduced linear system */
1609566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&gpcg->R));
1619566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&gpcg->DXFree));
1629566063dSJacob Faibussowitsch       PetscCall(TaoVecGetSubVec(tao->gradient, gpcg->Free_Local, tao->subset_type, 0.0, &gpcg->R));
1639566063dSJacob Faibussowitsch       PetscCall(VecScale(gpcg->R, -1.0));
1649566063dSJacob Faibussowitsch       PetscCall(TaoVecGetSubVec(tao->stepdirection, gpcg->Free_Local, tao->subset_type, 0.0, &gpcg->DXFree));
1659566063dSJacob Faibussowitsch       PetscCall(VecSet(gpcg->DXFree, 0.0));
16674962610SAlp Dener 
1679566063dSJacob Faibussowitsch       PetscCall(TaoMatGetSubMat(tao->hessian, gpcg->Free_Local, gpcg->Work, tao->subset_type, &gpcg->Hsub));
16874962610SAlp Dener 
16974962610SAlp Dener       if (tao->hessian_pre == tao->hessian) {
1709566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&gpcg->Hsub_pre));
1719566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)gpcg->Hsub));
17274962610SAlp Dener         gpcg->Hsub_pre = gpcg->Hsub;
17374962610SAlp Dener       } else {
1749566063dSJacob Faibussowitsch         PetscCall(TaoMatGetSubMat(tao->hessian, gpcg->Free_Local, gpcg->Work, tao->subset_type, &gpcg->Hsub_pre));
17574962610SAlp Dener       }
17674962610SAlp Dener 
1779566063dSJacob Faibussowitsch       PetscCall(KSPReset(tao->ksp));
1789566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(tao->ksp, gpcg->Hsub, gpcg->Hsub_pre));
17974962610SAlp Dener 
1809566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, gpcg->R, gpcg->DXFree));
1819566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &its));
18274962610SAlp Dener       tao->ksp_its += its;
18374962610SAlp Dener       tao->ksp_tot_its += its;
1849566063dSJacob Faibussowitsch       PetscCall(VecSet(tao->stepdirection, 0.0));
1859566063dSJacob Faibussowitsch       PetscCall(VecISAXPY(tao->stepdirection, gpcg->Free_Local, 1.0, gpcg->DXFree));
18674962610SAlp Dener 
1879566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
1889566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
18974962610SAlp Dener       f_new = f;
1909566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, &stepsize, &ls_status));
19174962610SAlp Dener 
19274962610SAlp Dener       actred = f_new - f;
19374962610SAlp Dener 
19474962610SAlp Dener       /* Evaluate the function and gradient at the new point */
1959566063dSJacob Faibussowitsch       PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->PG));
1969566063dSJacob Faibussowitsch       PetscCall(VecNorm(gpcg->PG, NORM_2, &gnorm));
19774962610SAlp Dener       f = f_new;
1989566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&gpcg->Free_Local));
1999566063dSJacob Faibussowitsch       PetscCall(VecWhichInactive(tao->XL, tao->solution, tao->gradient, tao->XU, PETSC_TRUE, &gpcg->Free_Local));
20074962610SAlp Dener     } else {
2019371c9d4SSatish Balay       actred     = 0;
2029371c9d4SSatish Balay       gpcg->step = 1.0;
20374962610SAlp Dener       /* if there were no free variables, no cg method */
20474962610SAlp Dener     }
20574962610SAlp Dener 
20674962610SAlp Dener     tao->niter++;
2079371c9d4SSatish Balay     gpcg->f      = f;
2089371c9d4SSatish Balay     gpcg->gnorm  = gnorm;
2099371c9d4SSatish Balay     gpcg->actred = actred;
2109566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, gpcg->gnorm, 0.0, tao->ksp_its));
2119566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, gpcg->gnorm, 0.0, tao->step));
212dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
21374962610SAlp Dener     if (tao->reason != TAO_CONTINUE_ITERATING) break;
21474962610SAlp Dener   } /* END MAIN LOOP  */
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21674962610SAlp Dener }
21774962610SAlp Dener 
GPCGGradProjections(Tao tao)218d71ae5a4SJacob Faibussowitsch static PetscErrorCode GPCGGradProjections(Tao tao)
219d71ae5a4SJacob Faibussowitsch {
22074962610SAlp Dener   TAO_GPCG                    *gpcg = (TAO_GPCG *)tao->data;
22174962610SAlp Dener   PetscInt                     i;
22274962610SAlp Dener   PetscReal                    actred = -1.0, actred_max = 0.0, gAg, gtg = gpcg->gnorm, alpha;
22374962610SAlp Dener   PetscReal                    f_new, gdx, stepsize;
22474962610SAlp Dener   Vec                          DX = tao->stepdirection, XL = tao->XL, XU = tao->XU, Work = gpcg->Work;
22574962610SAlp Dener   Vec                          X = tao->solution, G = tao->gradient;
22674962610SAlp Dener   TaoLineSearchConvergedReason lsflag = TAOLINESEARCH_CONTINUE_ITERATING;
22774962610SAlp Dener 
22874962610SAlp Dener   /*
22974962610SAlp Dener      The free, active, and binding variables should be already identified
23074962610SAlp Dener   */
23174962610SAlp Dener   PetscFunctionBegin;
23274962610SAlp Dener   for (i = 0; i < gpcg->maxgpits; i++) {
23374962610SAlp Dener     if (-actred <= (gpcg->pg_ftol) * actred_max) break;
2349566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(G, X, XL, XU, DX));
2359566063dSJacob Faibussowitsch     PetscCall(VecScale(DX, -1.0));
2369566063dSJacob Faibussowitsch     PetscCall(VecDot(DX, G, &gdx));
23774962610SAlp Dener 
2389566063dSJacob Faibussowitsch     PetscCall(MatMult(tao->hessian, DX, Work));
2399566063dSJacob Faibussowitsch     PetscCall(VecDot(DX, Work, &gAg));
24074962610SAlp Dener 
24174962610SAlp Dener     gpcg->gp_iterates++;
24274962610SAlp Dener     gpcg->total_gp_its++;
24374962610SAlp Dener 
24474962610SAlp Dener     gtg = -gdx;
24574962610SAlp Dener     if (PetscAbsReal(gAg) == 0.0) {
24674962610SAlp Dener       alpha = 1.0;
24774962610SAlp Dener     } else {
24874962610SAlp Dener       alpha = PetscAbsReal(gtg / gAg);
24974962610SAlp Dener     }
2509566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, alpha));
25174962610SAlp Dener     f_new = gpcg->f;
2529566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, X, &f_new, G, DX, &stepsize, &lsflag));
25374962610SAlp Dener 
25474962610SAlp Dener     /* Update the iterate */
25574962610SAlp Dener     actred     = f_new - gpcg->f;
25674962610SAlp Dener     actred_max = PetscMax(actred_max, -(f_new - gpcg->f));
25774962610SAlp Dener     gpcg->f    = f_new;
2589566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&gpcg->Free_Local));
2599566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(XL, X, tao->gradient, XU, PETSC_TRUE, &gpcg->Free_Local));
26074962610SAlp Dener   }
26174962610SAlp Dener 
26274962610SAlp Dener   gpcg->gnorm = gtg;
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26474962610SAlp Dener } /* End gradient projections */
26574962610SAlp Dener 
TaoComputeDual_GPCG(Tao tao,Vec DXL,Vec DXU)266d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoComputeDual_GPCG(Tao tao, Vec DXL, Vec DXU)
267d71ae5a4SJacob Faibussowitsch {
26874962610SAlp Dener   TAO_GPCG *gpcg = (TAO_GPCG *)tao->data;
26974962610SAlp Dener 
27074962610SAlp Dener   PetscFunctionBegin;
2719566063dSJacob Faibussowitsch   PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->Work));
2729566063dSJacob Faibussowitsch   PetscCall(VecCopy(gpcg->Work, DXL));
2739566063dSJacob Faibussowitsch   PetscCall(VecAXPY(DXL, -1.0, tao->gradient));
2749566063dSJacob Faibussowitsch   PetscCall(VecSet(DXU, 0.0));
2759566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMax(DXL, DXL, DXU));
27674962610SAlp Dener 
2779566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, DXU));
2789566063dSJacob Faibussowitsch   PetscCall(VecAXPY(DXU, -1.0, gpcg->Work));
2799566063dSJacob Faibussowitsch   PetscCall(VecSet(gpcg->Work, 0.0));
2809566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMin(DXU, gpcg->Work, DXU));
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28274962610SAlp Dener }
28374962610SAlp Dener 
28474962610SAlp Dener /*MC
28574962610SAlp Dener   TAOGPCG - gradient projected conjugate gradient algorithm is an active-set
28674962610SAlp Dener         conjugate-gradient based method for bound-constrained minimization
28774962610SAlp Dener 
28874962610SAlp Dener   Options Database Keys:
28974962610SAlp Dener + -tao_gpcg_maxpgits - maximum number of gradient projections for GPCG iterate
29074962610SAlp Dener - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets
29174962610SAlp Dener 
29274962610SAlp Dener   Level: beginner
29374962610SAlp Dener M*/
TaoCreate_GPCG(Tao tao)294d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_GPCG(Tao tao)
295d71ae5a4SJacob Faibussowitsch {
29674962610SAlp Dener   TAO_GPCG *gpcg;
29774962610SAlp Dener 
29874962610SAlp Dener   PetscFunctionBegin;
29974962610SAlp Dener   tao->ops->setup          = TaoSetup_GPCG;
30074962610SAlp Dener   tao->ops->solve          = TaoSolve_GPCG;
30174962610SAlp Dener   tao->ops->view           = TaoView_GPCG;
30274962610SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_GPCG;
30374962610SAlp Dener   tao->ops->destroy        = TaoDestroy_GPCG;
30474962610SAlp Dener   tao->ops->computedual    = TaoComputeDual_GPCG;
30574962610SAlp Dener 
3064dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&gpcg));
30774962610SAlp Dener   tao->data = (void *)gpcg;
30874962610SAlp Dener 
30974962610SAlp Dener   /* Override default settings (unless already changed) */
310606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
311606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 500);
312606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_funcs, 100000);
313606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1e-6 : 1e-12);
314606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, grtol, PetscDefined(USE_REAL_SINGLE) ? 1e-6 : 1e-12);
31574962610SAlp Dener 
31674962610SAlp Dener   /* Initialize pointers and variables */
31774962610SAlp Dener   gpcg->n        = 0;
31874962610SAlp Dener   gpcg->maxgpits = 8;
31974962610SAlp Dener   gpcg->pg_ftol  = 0.1;
32074962610SAlp Dener 
32174962610SAlp Dener   gpcg->gp_iterates  = 0; /* Cumulative number */
32274962610SAlp Dener   gpcg->total_gp_its = 0;
32374962610SAlp Dener 
32474962610SAlp Dener   /* Initialize pointers and variables */
32574962610SAlp Dener   gpcg->n_bind      = 0;
32674962610SAlp Dener   gpcg->n_free      = 0;
32774962610SAlp Dener   gpcg->n_upper     = 0;
32874962610SAlp Dener   gpcg->n_lower     = 0;
32974962610SAlp Dener   gpcg->subset_type = TAO_SUBSET_MASK;
33074962610SAlp Dener   gpcg->Hsub        = NULL;
33174962610SAlp Dener   gpcg->Hsub_pre    = NULL;
33274962610SAlp Dener 
3339566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
3359566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
3369566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPNASH));
33774962610SAlp Dener 
3389566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3399566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3409566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHGPCG));
3419566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, GPCGObjectiveAndGradient, tao));
3429566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34474962610SAlp Dener }
345