1aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/bmrm/bmrm.h>
2a7e14dcfSSatish Balay
3a7e14dcfSSatish Balay static PetscErrorCode init_df_solver(TAO_DF *);
4a7e14dcfSSatish Balay static PetscErrorCode ensure_df_space(PetscInt, TAO_DF *);
5a7e14dcfSSatish Balay static PetscErrorCode destroy_df_solver(TAO_DF *);
60e660641SBarry Smith static PetscReal phi(PetscReal *, PetscInt, PetscReal, PetscReal *, PetscReal, PetscReal *, PetscReal *, PetscReal *);
70e660641SBarry Smith static PetscInt project(PetscInt, PetscReal *, PetscReal, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, TAO_DF *);
8a7e14dcfSSatish Balay
solve(TAO_DF * df)9ad349883SBarry Smith static PetscErrorCode solve(TAO_DF *df)
10ad349883SBarry Smith {
11ad349883SBarry Smith PetscInt i, j, innerIter, it, it2, luv, info;
12ad349883SBarry Smith PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam = 0.0, lam_ext;
13ad349883SBarry Smith PetscReal DELTAsv, ProdDELTAsv;
14ad349883SBarry Smith PetscReal c, *tempQ;
15ad349883SBarry Smith PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol;
16ad349883SBarry Smith PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd;
17ad349883SBarry Smith PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk;
18ad349883SBarry Smith PetscReal **Q = df->Q, *f = df->f, *t = df->t;
19ad349883SBarry Smith PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv;
20ad349883SBarry Smith
21a748edf9SJed Brown PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected dim %" PetscInt_FMT " >= 0", dim);
22ad349883SBarry Smith /* variables for the adaptive nonmonotone linesearch */
23ad349883SBarry Smith PetscInt L, llast;
24ad349883SBarry Smith PetscReal fr, fbest, fv, fc, fv0;
25ad349883SBarry Smith
26ad349883SBarry Smith c = BMRM_INFTY;
27ad349883SBarry Smith
28ad349883SBarry Smith DELTAsv = EPS_SV;
29ad349883SBarry Smith if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F;
30ad349883SBarry Smith else ProdDELTAsv = EPS_SV;
31ad349883SBarry Smith
32ad349883SBarry Smith for (i = 0; i < dim; i++) tempv[i] = -x[i];
33ad349883SBarry Smith
34ad349883SBarry Smith lam_ext = 0.0;
35ad349883SBarry Smith
36ad349883SBarry Smith /* Project the initial solution */
37ad349883SBarry Smith project(dim, a, b, tempv, l, u, x, &lam_ext, df);
38ad349883SBarry Smith
39ad349883SBarry Smith /* Compute gradient
40ad349883SBarry Smith g = Q*x + f; */
41ad349883SBarry Smith
42ad349883SBarry Smith it = 0;
43ad349883SBarry Smith for (i = 0; i < dim; i++) {
44ad349883SBarry Smith if (PetscAbsReal(x[i]) > ProdDELTAsv) ipt[it++] = i;
45ad349883SBarry Smith }
46ad349883SBarry Smith
47ad349883SBarry Smith PetscCall(PetscArrayzero(t, dim));
48ad349883SBarry Smith for (i = 0; i < it; i++) {
49ad349883SBarry Smith tempQ = Q[ipt[i]];
50ad349883SBarry Smith for (j = 0; j < dim; j++) t[j] += (tempQ[j] * x[ipt[i]]);
51ad349883SBarry Smith }
52ad349883SBarry Smith for (i = 0; i < dim; i++) g[i] = t[i] + f[i];
53ad349883SBarry Smith
54ad349883SBarry Smith /* y = -(x_{k} - g_{k}) */
55ad349883SBarry Smith for (i = 0; i < dim; i++) y[i] = g[i] - x[i];
56ad349883SBarry Smith
57ad349883SBarry Smith /* Project x_{k} - g_{k} */
58ad349883SBarry Smith project(dim, a, b, y, l, u, tempv, &lam_ext, df);
59ad349883SBarry Smith
60ad349883SBarry Smith /* y = P(x_{k} - g_{k}) - x_{k} */
61ad349883SBarry Smith max = ALPHA_MIN;
62ad349883SBarry Smith for (i = 0; i < dim; i++) {
63ad349883SBarry Smith y[i] = tempv[i] - x[i];
64ad349883SBarry Smith if (PetscAbsReal(y[i]) > max) max = PetscAbsReal(y[i]);
65ad349883SBarry Smith }
66ad349883SBarry Smith
67ad349883SBarry Smith if (max < tol * 1e-3) return PETSC_SUCCESS;
68ad349883SBarry Smith
69ad349883SBarry Smith alpha = 1.0 / max;
70ad349883SBarry Smith
71ad349883SBarry Smith /* fv0 = f(x_{0}). Recall t = Q x_{k} */
72ad349883SBarry Smith fv0 = 0.0;
73ad349883SBarry Smith for (i = 0; i < dim; i++) fv0 += x[i] * (0.5 * t[i] + f[i]);
74ad349883SBarry Smith
75ad349883SBarry Smith /* adaptive nonmonotone linesearch */
76ad349883SBarry Smith L = 2;
77ad349883SBarry Smith fr = ALPHA_MAX;
78ad349883SBarry Smith fbest = fv0;
79ad349883SBarry Smith fc = fv0;
80ad349883SBarry Smith llast = 0;
81ad349883SBarry Smith akold = bkold = 0.0;
82ad349883SBarry Smith
83ad349883SBarry Smith /* Iterator begins */
84ad349883SBarry Smith for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) {
85ad349883SBarry Smith /* tempv = -(x_{k} - alpha*g_{k}) */
86ad349883SBarry Smith for (i = 0; i < dim; i++) tempv[i] = alpha * g[i] - x[i];
87ad349883SBarry Smith
88ad349883SBarry Smith /* Project x_{k} - alpha*g_{k} */
89ad349883SBarry Smith project(dim, a, b, tempv, l, u, y, &lam_ext, df);
90ad349883SBarry Smith
91ad349883SBarry Smith /* gd = \inner{d_{k}}{g_{k}}
92ad349883SBarry Smith d = P(x_{k} - alpha*g_{k}) - x_{k}
93ad349883SBarry Smith */
94ad349883SBarry Smith gd = 0.0;
95ad349883SBarry Smith for (i = 0; i < dim; i++) {
96ad349883SBarry Smith d[i] = y[i] - x[i];
97ad349883SBarry Smith gd += d[i] * g[i];
98ad349883SBarry Smith }
99ad349883SBarry Smith
100ad349883SBarry Smith /* Gradient computation */
101ad349883SBarry Smith
102ad349883SBarry Smith /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */
103ad349883SBarry Smith
104ad349883SBarry Smith it = it2 = 0;
105ad349883SBarry Smith for (i = 0; i < dim; i++) {
106ad349883SBarry Smith if (PetscAbsReal(d[i]) > (ProdDELTAsv * 1.0e-2)) ipt[it++] = i;
107ad349883SBarry Smith }
108ad349883SBarry Smith for (i = 0; i < dim; i++) {
109ad349883SBarry Smith if (PetscAbsReal(y[i]) > ProdDELTAsv) ipt2[it2++] = i;
110ad349883SBarry Smith }
111ad349883SBarry Smith
112ad349883SBarry Smith PetscCall(PetscArrayzero(Qd, dim));
113ad349883SBarry Smith /* compute Qd = Q*d */
114ad349883SBarry Smith if (it < it2) {
115ad349883SBarry Smith for (i = 0; i < it; i++) {
116ad349883SBarry Smith tempQ = Q[ipt[i]];
117ad349883SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]);
118ad349883SBarry Smith }
119ad349883SBarry Smith } else { /* compute Qd = Q*y-t */
120ad349883SBarry Smith for (i = 0; i < it2; i++) {
121ad349883SBarry Smith tempQ = Q[ipt2[i]];
122ad349883SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]);
123ad349883SBarry Smith }
124ad349883SBarry Smith for (j = 0; j < dim; j++) Qd[j] -= t[j];
125ad349883SBarry Smith }
126ad349883SBarry Smith
127ad349883SBarry Smith /* ak = inner{d_{k}}{d_{k}} */
128ad349883SBarry Smith ak = 0.0;
129ad349883SBarry Smith for (i = 0; i < dim; i++) ak += d[i] * d[i];
130ad349883SBarry Smith
131ad349883SBarry Smith bk = 0.0;
132ad349883SBarry Smith for (i = 0; i < dim; i++) bk += d[i] * Qd[i];
133ad349883SBarry Smith
134ad349883SBarry Smith if (bk > EPS * ak && gd < 0.0) lamnew = -gd / bk;
135ad349883SBarry Smith else lamnew = 1.0;
136ad349883SBarry Smith
137ad349883SBarry Smith /* fv is computing f(x_{k} + d_{k}) */
138ad349883SBarry Smith fv = 0.0;
139ad349883SBarry Smith for (i = 0; i < dim; i++) {
140ad349883SBarry Smith xplus[i] = x[i] + d[i];
141ad349883SBarry Smith tplus[i] = t[i] + Qd[i];
142ad349883SBarry Smith fv += xplus[i] * (0.5 * tplus[i] + f[i]);
143ad349883SBarry Smith }
144ad349883SBarry Smith
145ad349883SBarry Smith /* fr is fref */
146ad349883SBarry Smith if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)) {
147ad349883SBarry Smith fv = 0.0;
148ad349883SBarry Smith for (i = 0; i < dim; i++) {
149ad349883SBarry Smith xplus[i] = x[i] + lamnew * d[i];
150ad349883SBarry Smith tplus[i] = t[i] + lamnew * Qd[i];
151ad349883SBarry Smith fv += xplus[i] * (0.5 * tplus[i] + f[i]);
152ad349883SBarry Smith }
153ad349883SBarry Smith }
154ad349883SBarry Smith
155ad349883SBarry Smith for (i = 0; i < dim; i++) {
156ad349883SBarry Smith sk[i] = xplus[i] - x[i];
157ad349883SBarry Smith yk[i] = tplus[i] - t[i];
158ad349883SBarry Smith x[i] = xplus[i];
159ad349883SBarry Smith t[i] = tplus[i];
160ad349883SBarry Smith g[i] = t[i] + f[i];
161ad349883SBarry Smith }
162ad349883SBarry Smith
163ad349883SBarry Smith /* update the line search control parameters */
164ad349883SBarry Smith if (fv < fbest) {
165ad349883SBarry Smith fbest = fv;
166ad349883SBarry Smith fc = fv;
167ad349883SBarry Smith llast = 0;
168ad349883SBarry Smith } else {
169ad349883SBarry Smith fc = (fc > fv ? fc : fv);
170ad349883SBarry Smith llast++;
171ad349883SBarry Smith if (llast == L) {
172ad349883SBarry Smith fr = fc;
173ad349883SBarry Smith fc = fv;
174ad349883SBarry Smith llast = 0;
175ad349883SBarry Smith }
176ad349883SBarry Smith }
177ad349883SBarry Smith
178ad349883SBarry Smith ak = bk = 0.0;
179ad349883SBarry Smith for (i = 0; i < dim; i++) {
180ad349883SBarry Smith ak += sk[i] * sk[i];
181ad349883SBarry Smith bk += sk[i] * yk[i];
182ad349883SBarry Smith }
183ad349883SBarry Smith
184ad349883SBarry Smith if (bk <= EPS * ak) alpha = ALPHA_MAX;
185ad349883SBarry Smith else {
186ad349883SBarry Smith if (bkold < EPS * akold) alpha = ak / bk;
187ad349883SBarry Smith else alpha = (akold + ak) / (bkold + bk);
188ad349883SBarry Smith
189ad349883SBarry Smith if (alpha > ALPHA_MAX) alpha = ALPHA_MAX;
190ad349883SBarry Smith else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN;
191ad349883SBarry Smith }
192ad349883SBarry Smith
193ad349883SBarry Smith akold = ak;
194ad349883SBarry Smith bkold = bk;
195ad349883SBarry Smith
196ad349883SBarry Smith /* stopping criterion based on KKT conditions */
197ad349883SBarry Smith /* at optimal, gradient of lagrangian w.r.t. x is zero */
198ad349883SBarry Smith
199ad349883SBarry Smith bk = 0.0;
200ad349883SBarry Smith for (i = 0; i < dim; i++) bk += x[i] * x[i];
201ad349883SBarry Smith
202ad349883SBarry Smith if (PetscSqrtReal(ak) < tol * 10 * PetscSqrtReal(bk)) {
203ad349883SBarry Smith it = 0;
204ad349883SBarry Smith luv = 0;
205ad349883SBarry Smith kktlam = 0.0;
206ad349883SBarry Smith for (i = 0; i < dim; i++) {
207ad349883SBarry Smith /* x[i] is active hence lagrange multipliers for box constraints
208ad349883SBarry Smith are zero. The lagrange multiplier for ineq. const. is then
209ad349883SBarry Smith defined as below
210ad349883SBarry Smith */
211ad349883SBarry Smith if ((x[i] > DELTAsv) && (x[i] < c - DELTAsv)) {
212ad349883SBarry Smith ipt[it++] = i;
213ad349883SBarry Smith kktlam = kktlam - a[i] * g[i];
214ad349883SBarry Smith } else uv[luv++] = i;
215ad349883SBarry Smith }
216ad349883SBarry Smith
217ad349883SBarry Smith if (it == 0 && PetscSqrtReal(ak) < tol * 0.5 * PetscSqrtReal(bk)) return PETSC_SUCCESS;
218ad349883SBarry Smith else {
219ad349883SBarry Smith kktlam = kktlam / it;
220ad349883SBarry Smith info = 1;
221ad349883SBarry Smith for (i = 0; i < it; i++) {
222ad349883SBarry Smith if (PetscAbsReal(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) {
223ad349883SBarry Smith info = 0;
224ad349883SBarry Smith break;
225ad349883SBarry Smith }
226ad349883SBarry Smith }
227ad349883SBarry Smith if (info == 1) {
228ad349883SBarry Smith for (i = 0; i < luv; i++) {
229ad349883SBarry Smith if (x[uv[i]] <= DELTAsv) {
230ad349883SBarry Smith /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may
231ad349883SBarry Smith not be zero. So, the gradient without beta is > 0
232ad349883SBarry Smith */
233ad349883SBarry Smith if (g[uv[i]] + kktlam * a[uv[i]] < -tol) {
234ad349883SBarry Smith info = 0;
235ad349883SBarry Smith break;
236ad349883SBarry Smith }
237ad349883SBarry Smith } else {
238ad349883SBarry Smith /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may
239ad349883SBarry Smith not be zero. So, the gradient without eta is < 0
240ad349883SBarry Smith */
241ad349883SBarry Smith if (g[uv[i]] + kktlam * a[uv[i]] > tol) {
242ad349883SBarry Smith info = 0;
243ad349883SBarry Smith break;
244ad349883SBarry Smith }
245ad349883SBarry Smith }
246ad349883SBarry Smith }
247ad349883SBarry Smith }
248ad349883SBarry Smith
249ad349883SBarry Smith if (info == 1) return PETSC_SUCCESS;
250ad349883SBarry Smith }
251ad349883SBarry Smith }
252ad349883SBarry Smith }
253ad349883SBarry Smith return PETSC_SUCCESS;
254ad349883SBarry Smith }
255ad349883SBarry Smith
256a7e14dcfSSatish Balay /* The main solver function
257a7e14dcfSSatish Balay
258a7e14dcfSSatish Balay f = Remp(W) This is what the user provides us from the application layer
259a7e14dcfSSatish Balay So the ComputeGradient function for instance should get us back the subgradient of Remp(W)
260a7e14dcfSSatish Balay
261a7e14dcfSSatish Balay Regularizer assumed to be L2 norm = lambda*0.5*W'W ()
262a7e14dcfSSatish Balay */
263a7e14dcfSSatish Balay
make_grad_node(Vec X,Vec_Chain ** p)264d71ae5a4SJacob Faibussowitsch static PetscErrorCode make_grad_node(Vec X, Vec_Chain **p)
265d71ae5a4SJacob Faibussowitsch {
266a7e14dcfSSatish Balay PetscFunctionBegin;
2679566063dSJacob Faibussowitsch PetscCall(PetscNew(p));
2689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &(*p)->V));
2699566063dSJacob Faibussowitsch PetscCall(VecCopy(X, (*p)->V));
2706c23d075SBarry Smith (*p)->next = NULL;
2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
272a7e14dcfSSatish Balay }
273a7e14dcfSSatish Balay
destroy_grad_list(Vec_Chain * head)274d71ae5a4SJacob Faibussowitsch static PetscErrorCode destroy_grad_list(Vec_Chain *head)
275d71ae5a4SJacob Faibussowitsch {
276a7e14dcfSSatish Balay Vec_Chain *p = head->next, *q;
277a7e14dcfSSatish Balay
278a7e14dcfSSatish Balay PetscFunctionBegin;
279a7e14dcfSSatish Balay while (p) {
280a7e14dcfSSatish Balay q = p->next;
2819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p->V));
2829566063dSJacob Faibussowitsch PetscCall(PetscFree(p));
283a7e14dcfSSatish Balay p = q;
284a7e14dcfSSatish Balay }
2856c23d075SBarry Smith head->next = NULL;
2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
287a7e14dcfSSatish Balay }
288a7e14dcfSSatish Balay
TaoSolve_BMRM(Tao tao)289d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BMRM(Tao tao)
290d71ae5a4SJacob Faibussowitsch {
291a7e14dcfSSatish Balay TAO_DF df;
292a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM *)tao->data;
293a7e14dcfSSatish Balay
294a7e14dcfSSatish Balay /* Values and pointers to parts of the optimization problem */
295a7e14dcfSSatish Balay PetscReal f = 0.0;
296a7e14dcfSSatish Balay Vec W = tao->solution;
297a7e14dcfSSatish Balay Vec G = tao->gradient;
298a7e14dcfSSatish Balay PetscReal lambda;
299a7e14dcfSSatish Balay PetscReal bt;
300a7e14dcfSSatish Balay Vec_Chain grad_list, *tail_glist, *pgrad;
301a7e14dcfSSatish Balay PetscInt i;
302a7e14dcfSSatish Balay PetscMPIInt rank;
303a7e14dcfSSatish Balay
304e4cb33bbSBarry Smith /* Used in converged criteria check */
305a7e14dcfSSatish Balay PetscReal reg;
3067fb8a5e4SKarl Rupp PetscReal jtwt = 0.0, max_jtwt, pre_epsilon, epsilon, jw, min_jw;
307a7e14dcfSSatish Balay PetscReal innerSolverTol;
308ba4b436cSBarry Smith MPI_Comm comm;
309a7e14dcfSSatish Balay
310a7e14dcfSSatish Balay PetscFunctionBegin;
3119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
3129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
313a7e14dcfSSatish Balay lambda = bmrm->lambda;
314a7e14dcfSSatish Balay
315a7e14dcfSSatish Balay /* Check Stopping Condition */
316a7e14dcfSSatish Balay tao->step = 1.0;
317a7e14dcfSSatish Balay max_jtwt = -BMRM_INFTY;
318a7e14dcfSSatish Balay min_jw = BMRM_INFTY;
319a7e14dcfSSatish Balay innerSolverTol = 1.0;
320a7e14dcfSSatish Balay epsilon = 0.0;
321a7e14dcfSSatish Balay
322dd400576SPatrick Sanan if (rank == 0) {
3239566063dSJacob Faibussowitsch PetscCall(init_df_solver(&df));
324a7e14dcfSSatish Balay grad_list.next = NULL;
325a7e14dcfSSatish Balay tail_glist = &grad_list;
326a7e14dcfSSatish Balay }
327a7e14dcfSSatish Balay
328a7e14dcfSSatish Balay df.tol = 1e-6;
3293ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
330a7e14dcfSSatish Balay
331a7e14dcfSSatish Balay /*-----------------Algorithm Begins------------------------*/
332a7e14dcfSSatish Balay /* make the scatter */
3339566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w));
3349566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bmrm->local_w));
3359566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bmrm->local_w));
336a7e14dcfSSatish Balay
337a7e14dcfSSatish Balay /* NOTE: In application pass the sub-gradient of Remp(W) */
3389566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, W, &f, G));
3399566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, 1.0, 0.0, tao->ksp_its));
3409566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, 1.0, 0.0, tao->step));
341dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
3423ecd9318SAlp Dener
3433ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) {
344e1e80dc8SAlp Dener /* Call general purpose update function */
345dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
346e1e80dc8SAlp Dener
347a7e14dcfSSatish Balay /* compute bt = Remp(Wt-1) - <Wt-1, At> */
3489566063dSJacob Faibussowitsch PetscCall(VecDot(W, G, &bt));
349a7e14dcfSSatish Balay bt = f - bt;
350a7e14dcfSSatish Balay
3519dddd249SSatish Balay /* First gather the gradient to the rank-0 node */
3529566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD));
3539566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD));
354a7e14dcfSSatish Balay
355a7e14dcfSSatish Balay /* Bring up the inner solver */
356dd400576SPatrick Sanan if (rank == 0) {
3579566063dSJacob Faibussowitsch PetscCall(ensure_df_space(tao->niter + 1, &df));
3589566063dSJacob Faibussowitsch PetscCall(make_grad_node(bmrm->local_w, &pgrad));
359a7e14dcfSSatish Balay tail_glist->next = pgrad;
360a7e14dcfSSatish Balay tail_glist = pgrad;
361a7e14dcfSSatish Balay
3628931d482SJason Sarich df.a[tao->niter] = 1.0;
3638931d482SJason Sarich df.f[tao->niter] = -bt;
3648931d482SJason Sarich df.u[tao->niter] = 1.0;
3658931d482SJason Sarich df.l[tao->niter] = 0.0;
366a7e14dcfSSatish Balay
367a7e14dcfSSatish Balay /* set up the Q */
368a7e14dcfSSatish Balay pgrad = grad_list.next;
3698931d482SJason Sarich for (i = 0; i <= tao->niter; i++) {
3703c859ba3SBarry Smith PetscCheck(pgrad, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Assert that there are at least tao->niter+1 pgrad available");
3719566063dSJacob Faibussowitsch PetscCall(VecDot(pgrad->V, bmrm->local_w, ®));
3728931d482SJason Sarich df.Q[i][tao->niter] = df.Q[tao->niter][i] = reg / lambda;
373a7e14dcfSSatish Balay pgrad = pgrad->next;
374a7e14dcfSSatish Balay }
375a7e14dcfSSatish Balay
3768931d482SJason Sarich if (tao->niter > 0) {
3778931d482SJason Sarich df.x[tao->niter] = 0.0;
3789566063dSJacob Faibussowitsch PetscCall(solve(&df));
3799371c9d4SSatish Balay } else df.x[0] = 1.0;
380a7e14dcfSSatish Balay
381a7e14dcfSSatish Balay /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */
382a7e14dcfSSatish Balay jtwt = 0.0;
3839566063dSJacob Faibussowitsch PetscCall(VecSet(bmrm->local_w, 0.0));
384a7e14dcfSSatish Balay pgrad = grad_list.next;
3858931d482SJason Sarich for (i = 0; i <= tao->niter; i++) {
386a7e14dcfSSatish Balay jtwt -= df.x[i] * df.f[i];
3879566063dSJacob Faibussowitsch PetscCall(VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V));
388a7e14dcfSSatish Balay pgrad = pgrad->next;
389a7e14dcfSSatish Balay }
390a7e14dcfSSatish Balay
3919566063dSJacob Faibussowitsch PetscCall(VecNorm(bmrm->local_w, NORM_2, ®));
392a7e14dcfSSatish Balay reg = 0.5 * lambda * reg * reg;
393a7e14dcfSSatish Balay jtwt -= reg;
394a7e14dcfSSatish Balay } /* end if rank == 0 */
395a7e14dcfSSatish Balay
396a7e14dcfSSatish Balay /* scatter the new W to all nodes */
3979566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(bmrm->scatter, bmrm->local_w, W, INSERT_VALUES, SCATTER_REVERSE));
3989566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(bmrm->scatter, bmrm->local_w, W, INSERT_VALUES, SCATTER_REVERSE));
399a7e14dcfSSatish Balay
4009566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, W, &f, G));
401a7e14dcfSSatish Balay
4029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&jtwt, 1, MPIU_REAL, 0, comm));
4039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(®, 1, MPIU_REAL, 0, comm));
404a7e14dcfSSatish Balay
405a7e14dcfSSatish Balay jw = reg + f; /* J(w) = regularizer + Remp(w) */
4060e660641SBarry Smith if (jw < min_jw) min_jw = jw;
4070e660641SBarry Smith if (jtwt > max_jtwt) max_jtwt = jtwt;
408a7e14dcfSSatish Balay
409a7e14dcfSSatish Balay pre_epsilon = epsilon;
410a7e14dcfSSatish Balay epsilon = min_jw - jtwt;
411a7e14dcfSSatish Balay
412dd400576SPatrick Sanan if (rank == 0) {
4130e660641SBarry Smith if (innerSolverTol > epsilon) innerSolverTol = epsilon;
4140e660641SBarry Smith else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7;
415a7e14dcfSSatish Balay
416a7e14dcfSSatish Balay /* if the annealing doesn't work well, lower the inner solver tolerance */
4170e660641SBarry Smith if (pre_epsilon < epsilon) innerSolverTol *= 0.2;
418a7e14dcfSSatish Balay
419a7e14dcfSSatish Balay df.tol = innerSolverTol * 0.5;
420a7e14dcfSSatish Balay }
421a7e14dcfSSatish Balay
4228931d482SJason Sarich tao->niter++;
4239566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, min_jw, epsilon, 0.0, tao->ksp_its));
4249566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, min_jw, epsilon, 0.0, tao->step));
425dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
426a7e14dcfSSatish Balay }
427a7e14dcfSSatish Balay
428a7e14dcfSSatish Balay /* free all the memory */
429dd400576SPatrick Sanan if (rank == 0) {
4309566063dSJacob Faibussowitsch PetscCall(destroy_grad_list(&grad_list));
4319566063dSJacob Faibussowitsch PetscCall(destroy_df_solver(&df));
432a7e14dcfSSatish Balay }
433a7e14dcfSSatish Balay
4349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmrm->local_w));
4359566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&bmrm->scatter));
4363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
437a7e14dcfSSatish Balay }
438a7e14dcfSSatish Balay
TaoSetup_BMRM(Tao tao)439d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_BMRM(Tao tao)
440d71ae5a4SJacob Faibussowitsch {
441a7e14dcfSSatish Balay PetscFunctionBegin;
442a7e14dcfSSatish Balay /* Allocate some arrays */
4439566063dSJacob Faibussowitsch if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
445a7e14dcfSSatish Balay }
446a7e14dcfSSatish Balay
TaoDestroy_BMRM(Tao tao)447d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BMRM(Tao tao)
448d71ae5a4SJacob Faibussowitsch {
449a7e14dcfSSatish Balay PetscFunctionBegin;
4509566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
4513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
452a7e14dcfSSatish Balay }
453a7e14dcfSSatish Balay
TaoSetFromOptions_BMRM(Tao tao,PetscOptionItems PetscOptionsObject)454*ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BMRM(Tao tao, PetscOptionItems PetscOptionsObject)
455d71ae5a4SJacob Faibussowitsch {
456a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM *)tao->data;
457a7e14dcfSSatish Balay
458a7e14dcfSSatish Balay PetscFunctionBegin;
459d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "BMRM for regularized risk minimization");
4609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight", "", 100, &bmrm->lambda, NULL));
461d0609cedSBarry Smith PetscOptionsHeadEnd();
4623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
463a7e14dcfSSatish Balay }
464a7e14dcfSSatish Balay
TaoView_BMRM(Tao tao,PetscViewer viewer)465d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer)
466d71ae5a4SJacob Faibussowitsch {
467a7e14dcfSSatish Balay PetscBool isascii;
468a7e14dcfSSatish Balay
469a7e14dcfSSatish Balay PetscFunctionBegin;
4709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
471a7e14dcfSSatish Balay if (isascii) {
4729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
4739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
474a7e14dcfSSatish Balay }
4753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
476a7e14dcfSSatish Balay }
477a7e14dcfSSatish Balay
4781522df2eSJason Sarich /*MC
4791522df2eSJason Sarich TAOBMRM - bundle method for regularized risk minimization
4801522df2eSJason Sarich
4811522df2eSJason Sarich Options Database Keys:
4821522df2eSJason Sarich . - tao_bmrm_lambda - regulariser weight
4831522df2eSJason Sarich
4841eb8069cSJason Sarich Level: beginner
4851522df2eSJason Sarich M*/
4861522df2eSJason Sarich
TaoCreate_BMRM(Tao tao)487d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BMRM(Tao tao)
488d71ae5a4SJacob Faibussowitsch {
489a7e14dcfSSatish Balay TAO_BMRM *bmrm;
490a7e14dcfSSatish Balay
491a7e14dcfSSatish Balay PetscFunctionBegin;
492a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BMRM;
493a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BMRM;
494a7e14dcfSSatish Balay tao->ops->view = TaoView_BMRM;
495a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BMRM;
496a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BMRM;
497a7e14dcfSSatish Balay
4984dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bmrm));
499a7e14dcfSSatish Balay bmrm->lambda = 1.0;
500a7e14dcfSSatish Balay tao->data = (void *)bmrm;
501a7e14dcfSSatish Balay
5026552cf8aSJason Sarich /* Override default settings (unless already changed) */
503606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
504606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 2000);
505606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 4000);
506606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, gatol, 1.0e-12);
507606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, grtol, 1.0e-12);
5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
509a7e14dcfSSatish Balay }
510a7e14dcfSSatish Balay
init_df_solver(TAO_DF * df)51166976f2fSJacob Faibussowitsch static PetscErrorCode init_df_solver(TAO_DF *df)
512d71ae5a4SJacob Faibussowitsch {
513a7e14dcfSSatish Balay PetscInt i, n = INCRE_DIM;
514a7e14dcfSSatish Balay
515a7e14dcfSSatish Balay PetscFunctionBegin;
516a7e14dcfSSatish Balay /* default values */
517a7e14dcfSSatish Balay df->maxProjIter = 200;
518a7e14dcfSSatish Balay df->maxPGMIter = 300000;
519a7e14dcfSSatish Balay df->b = 1.0;
520a7e14dcfSSatish Balay
521a7e14dcfSSatish Balay /* memory space required by Dai-Fletcher */
522a7e14dcfSSatish Balay df->cur_num_cp = n;
5239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->f));
5249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->a));
5259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->l));
5269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->u));
5279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->x));
5289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Q));
529a7e14dcfSSatish Balay
53048a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscMalloc1(n, &df->Q[i]));
531a7e14dcfSSatish Balay
5329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->g));
5339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->y));
5349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tempv));
5359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->d));
5369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Qd));
5379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->t));
5389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->xplus));
5399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tplus));
5409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->sk));
5419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->yk));
542a7e14dcfSSatish Balay
5439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt));
5449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt2));
5459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->uv));
5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
547a7e14dcfSSatish Balay }
548a7e14dcfSSatish Balay
ensure_df_space(PetscInt dim,TAO_DF * df)54966976f2fSJacob Faibussowitsch static PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df)
550d71ae5a4SJacob Faibussowitsch {
551a7e14dcfSSatish Balay PetscReal *tmp, **tmp_Q;
552a7e14dcfSSatish Balay PetscInt i, n, old_n;
553a7e14dcfSSatish Balay
554a7e14dcfSSatish Balay PetscFunctionBegin;
55553506e15SBarry Smith df->dim = dim;
5563ba16761SJacob Faibussowitsch if (dim <= df->cur_num_cp) PetscFunctionReturn(PETSC_SUCCESS);
557a7e14dcfSSatish Balay
558a7e14dcfSSatish Balay old_n = df->cur_num_cp;
559a7e14dcfSSatish Balay df->cur_num_cp += INCRE_DIM;
560a7e14dcfSSatish Balay n = df->cur_num_cp;
561a7e14dcfSSatish Balay
562a7e14dcfSSatish Balay /* memory space required by dai-fletcher */
5639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp));
5649566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->f, old_n));
5659566063dSJacob Faibussowitsch PetscCall(PetscFree(df->f));
566a7e14dcfSSatish Balay df->f = tmp;
567a7e14dcfSSatish Balay
5689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp));
5699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->a, old_n));
5709566063dSJacob Faibussowitsch PetscCall(PetscFree(df->a));
571a7e14dcfSSatish Balay df->a = tmp;
572a7e14dcfSSatish Balay
5739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp));
5749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->l, old_n));
5759566063dSJacob Faibussowitsch PetscCall(PetscFree(df->l));
576a7e14dcfSSatish Balay df->l = tmp;
577a7e14dcfSSatish Balay
5789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp));
5799566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->u, old_n));
5809566063dSJacob Faibussowitsch PetscCall(PetscFree(df->u));
581a7e14dcfSSatish Balay df->u = tmp;
582a7e14dcfSSatish Balay
5839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp));
5849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->x, old_n));
5859566063dSJacob Faibussowitsch PetscCall(PetscFree(df->x));
586a7e14dcfSSatish Balay df->x = tmp;
587a7e14dcfSSatish Balay
5889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp_Q));
58953506e15SBarry Smith for (i = 0; i < n; i++) {
5909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp_Q[i]));
59153506e15SBarry Smith if (i < old_n) {
5929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp_Q[i], df->Q[i], old_n));
5939566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q[i]));
594a7e14dcfSSatish Balay }
595a7e14dcfSSatish Balay }
596a7e14dcfSSatish Balay
5979566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q));
598a7e14dcfSSatish Balay df->Q = tmp_Q;
599a7e14dcfSSatish Balay
6009566063dSJacob Faibussowitsch PetscCall(PetscFree(df->g));
6019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->g));
602a7e14dcfSSatish Balay
6039566063dSJacob Faibussowitsch PetscCall(PetscFree(df->y));
6049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->y));
605a7e14dcfSSatish Balay
6069566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tempv));
6079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tempv));
608a7e14dcfSSatish Balay
6099566063dSJacob Faibussowitsch PetscCall(PetscFree(df->d));
6109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->d));
611a7e14dcfSSatish Balay
6129566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Qd));
6139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Qd));
614a7e14dcfSSatish Balay
6159566063dSJacob Faibussowitsch PetscCall(PetscFree(df->t));
6169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->t));
617a7e14dcfSSatish Balay
6189566063dSJacob Faibussowitsch PetscCall(PetscFree(df->xplus));
6199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->xplus));
620a7e14dcfSSatish Balay
6219566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tplus));
6229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tplus));
623a7e14dcfSSatish Balay
6249566063dSJacob Faibussowitsch PetscCall(PetscFree(df->sk));
6259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->sk));
626a7e14dcfSSatish Balay
6279566063dSJacob Faibussowitsch PetscCall(PetscFree(df->yk));
6289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->yk));
629a7e14dcfSSatish Balay
6309566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt));
6319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt));
632a7e14dcfSSatish Balay
6339566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt2));
6349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt2));
635a7e14dcfSSatish Balay
6369566063dSJacob Faibussowitsch PetscCall(PetscFree(df->uv));
6379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->uv));
6383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
639a7e14dcfSSatish Balay }
640a7e14dcfSSatish Balay
destroy_df_solver(TAO_DF * df)64166976f2fSJacob Faibussowitsch static PetscErrorCode destroy_df_solver(TAO_DF *df)
642d71ae5a4SJacob Faibussowitsch {
643a7e14dcfSSatish Balay PetscInt i;
6446c23d075SBarry Smith
645a7e14dcfSSatish Balay PetscFunctionBegin;
6469566063dSJacob Faibussowitsch PetscCall(PetscFree(df->f));
6479566063dSJacob Faibussowitsch PetscCall(PetscFree(df->a));
6489566063dSJacob Faibussowitsch PetscCall(PetscFree(df->l));
6499566063dSJacob Faibussowitsch PetscCall(PetscFree(df->u));
6509566063dSJacob Faibussowitsch PetscCall(PetscFree(df->x));
651a7e14dcfSSatish Balay
65248a46eb9SPierre Jolivet for (i = 0; i < df->cur_num_cp; i++) PetscCall(PetscFree(df->Q[i]));
6539566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q));
6549566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt));
6559566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt2));
6569566063dSJacob Faibussowitsch PetscCall(PetscFree(df->uv));
6579566063dSJacob Faibussowitsch PetscCall(PetscFree(df->g));
6589566063dSJacob Faibussowitsch PetscCall(PetscFree(df->y));
6599566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tempv));
6609566063dSJacob Faibussowitsch PetscCall(PetscFree(df->d));
6619566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Qd));
6629566063dSJacob Faibussowitsch PetscCall(PetscFree(df->t));
6639566063dSJacob Faibussowitsch PetscCall(PetscFree(df->xplus));
6649566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tplus));
6659566063dSJacob Faibussowitsch PetscCall(PetscFree(df->sk));
6669566063dSJacob Faibussowitsch PetscCall(PetscFree(df->yk));
6673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
668a7e14dcfSSatish Balay }
669a7e14dcfSSatish Balay
670a7e14dcfSSatish Balay /* Piecewise linear monotone target function for the Dai-Fletcher projector */
phi(PetscReal * x,PetscInt n,PetscReal lambda,PetscReal * a,PetscReal b,PetscReal * c,PetscReal * l,PetscReal * u)67166976f2fSJacob Faibussowitsch static PetscReal phi(PetscReal *x, PetscInt n, PetscReal lambda, PetscReal *a, PetscReal b, PetscReal *c, PetscReal *l, PetscReal *u)
672d71ae5a4SJacob Faibussowitsch {
673a7e14dcfSSatish Balay PetscReal r = 0.0;
674a7e14dcfSSatish Balay PetscInt i;
675a7e14dcfSSatish Balay
676a7e14dcfSSatish Balay for (i = 0; i < n; i++) {
677a7e14dcfSSatish Balay x[i] = -c[i] + lambda * a[i];
6786c23d075SBarry Smith if (x[i] > u[i]) x[i] = u[i];
6796c23d075SBarry Smith else if (x[i] < l[i]) x[i] = l[i];
680a7e14dcfSSatish Balay r += a[i] * x[i];
681a7e14dcfSSatish Balay }
682a7e14dcfSSatish Balay return r - b;
683a7e14dcfSSatish Balay }
684a7e14dcfSSatish Balay
685a7e14dcfSSatish Balay /** Modified Dai-Fletcher QP projector solves the problem:
686a7e14dcfSSatish Balay *
687a7e14dcfSSatish Balay * minimise 0.5*x'*x - c'*x
688a7e14dcfSSatish Balay * subj to a'*x = b
689a7e14dcfSSatish Balay * l \leq x \leq u
690a7e14dcfSSatish Balay *
691a7e14dcfSSatish Balay * \param c The point to be projected onto feasible set
692a7e14dcfSSatish Balay */
project(PetscInt n,PetscReal * a,PetscReal b,PetscReal * c,PetscReal * l,PetscReal * u,PetscReal * x,PetscReal * lam_ext,TAO_DF * df)69366976f2fSJacob Faibussowitsch static PetscInt project(PetscInt n, PetscReal *a, PetscReal b, PetscReal *c, PetscReal *l, PetscReal *u, PetscReal *x, PetscReal *lam_ext, TAO_DF *df)
694d71ae5a4SJacob Faibussowitsch {
695a7e14dcfSSatish Balay PetscReal lambda, lambdal, lambdau, dlambda, lambda_new;
696a7e14dcfSSatish Balay PetscReal r, rl, ru, s;
697a7e14dcfSSatish Balay PetscInt innerIter;
698a7e14dcfSSatish Balay PetscBool nonNegativeSlack = PETSC_FALSE;
699a7e14dcfSSatish Balay
700a7e14dcfSSatish Balay *lam_ext = 0;
701a7e14dcfSSatish Balay lambda = 0;
702a7e14dcfSSatish Balay dlambda = 0.5;
703a7e14dcfSSatish Balay innerIter = 1;
704a7e14dcfSSatish Balay
705a7e14dcfSSatish Balay /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b)
706a7e14dcfSSatish Balay *
707a7e14dcfSSatish Balay * Optimality conditions for \phi:
708a7e14dcfSSatish Balay *
709a7e14dcfSSatish Balay * 1. lambda <= 0
710a7e14dcfSSatish Balay * 2. r <= 0
711a7e14dcfSSatish Balay * 3. r*lambda == 0
712a7e14dcfSSatish Balay */
713a7e14dcfSSatish Balay
714a7e14dcfSSatish Balay /* Bracketing Phase */
715a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u);
716a7e14dcfSSatish Balay
7176c23d075SBarry Smith if (nonNegativeSlack) {
718a7e14dcfSSatish Balay /* inequality constraint, i.e., with \xi >= 0 constraint */
7193ba16761SJacob Faibussowitsch if (r < TOL_R) return PETSC_SUCCESS;
7206c23d075SBarry Smith } else {
721a7e14dcfSSatish Balay /* equality constraint ,i.e., without \xi >= 0 constraint */
7223ba16761SJacob Faibussowitsch if (PetscAbsReal(r) < TOL_R) return PETSC_SUCCESS;
723a7e14dcfSSatish Balay }
724a7e14dcfSSatish Balay
725a7e14dcfSSatish Balay if (r < 0.0) {
726a7e14dcfSSatish Balay lambdal = lambda;
727a7e14dcfSSatish Balay rl = r;
728a7e14dcfSSatish Balay lambda = lambda + dlambda;
729a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u);
730a7e14dcfSSatish Balay while (r < 0.0 && dlambda < BMRM_INFTY) {
731a7e14dcfSSatish Balay lambdal = lambda;
732a7e14dcfSSatish Balay s = rl / r - 1.0;
733a7e14dcfSSatish Balay if (s < 0.1) s = 0.1;
734a7e14dcfSSatish Balay dlambda = dlambda + dlambda / s;
735a7e14dcfSSatish Balay lambda = lambda + dlambda;
736a7e14dcfSSatish Balay rl = r;
737a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u);
738a7e14dcfSSatish Balay }
739a7e14dcfSSatish Balay lambdau = lambda;
740a7e14dcfSSatish Balay ru = r;
7416c23d075SBarry Smith } else {
742a7e14dcfSSatish Balay lambdau = lambda;
743a7e14dcfSSatish Balay ru = r;
744a7e14dcfSSatish Balay lambda = lambda - dlambda;
745a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u);
746a7e14dcfSSatish Balay while (r > 0.0 && dlambda > -BMRM_INFTY) {
747a7e14dcfSSatish Balay lambdau = lambda;
748a7e14dcfSSatish Balay s = ru / r - 1.0;
749a7e14dcfSSatish Balay if (s < 0.1) s = 0.1;
750a7e14dcfSSatish Balay dlambda = dlambda + dlambda / s;
751a7e14dcfSSatish Balay lambda = lambda - dlambda;
752a7e14dcfSSatish Balay ru = r;
753a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u);
754a7e14dcfSSatish Balay }
755a7e14dcfSSatish Balay lambdal = lambda;
756a7e14dcfSSatish Balay rl = r;
757a7e14dcfSSatish Balay }
758a7e14dcfSSatish Balay
7593387f462SBarry Smith PetscCheckAbort(PetscAbsReal(dlambda) <= BMRM_INFTY, PETSC_COMM_SELF, PETSC_ERR_PLIB, "L2N2_DaiFletcherPGM detected Infeasible QP problem!");
760a7e14dcfSSatish Balay
761ad540459SPierre Jolivet if (ru == 0) return innerIter;
762a7e14dcfSSatish Balay
763a7e14dcfSSatish Balay /* Secant Phase */
764a7e14dcfSSatish Balay s = 1.0 - rl / ru;
765a7e14dcfSSatish Balay dlambda = dlambda / s;
766a7e14dcfSSatish Balay lambda = lambdau - dlambda;
767a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u);
768a7e14dcfSSatish Balay
7699371c9d4SSatish Balay while (PetscAbsReal(r) > TOL_R && dlambda > TOL_LAM * (1.0 + PetscAbsReal(lambda)) && innerIter < df->maxProjIter) {
770a7e14dcfSSatish Balay innerIter++;
771a7e14dcfSSatish Balay if (r > 0.0) {
772a7e14dcfSSatish Balay if (s <= 2.0) {
773a7e14dcfSSatish Balay lambdau = lambda;
774a7e14dcfSSatish Balay ru = r;
775a7e14dcfSSatish Balay s = 1.0 - rl / ru;
776a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s;
777a7e14dcfSSatish Balay lambda = lambdau - dlambda;
77853506e15SBarry Smith } else {
779a7e14dcfSSatish Balay s = ru / r - 1.0;
780a7e14dcfSSatish Balay if (s < 0.1) s = 0.1;
781a7e14dcfSSatish Balay dlambda = (lambdau - lambda) / s;
782a7e14dcfSSatish Balay lambda_new = 0.75 * lambdal + 0.25 * lambda;
7839371c9d4SSatish Balay if (lambda_new < (lambda - dlambda)) lambda_new = lambda - dlambda;
784a7e14dcfSSatish Balay lambdau = lambda;
785a7e14dcfSSatish Balay ru = r;
786a7e14dcfSSatish Balay lambda = lambda_new;
787a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda);
788a7e14dcfSSatish Balay }
78953506e15SBarry Smith } else {
790a7e14dcfSSatish Balay if (s >= 2.0) {
791a7e14dcfSSatish Balay lambdal = lambda;
792a7e14dcfSSatish Balay rl = r;
793a7e14dcfSSatish Balay s = 1.0 - rl / ru;
794a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s;
795a7e14dcfSSatish Balay lambda = lambdau - dlambda;
79653506e15SBarry Smith } else {
797a7e14dcfSSatish Balay s = rl / r - 1.0;
798a7e14dcfSSatish Balay if (s < 0.1) s = 0.1;
799a7e14dcfSSatish Balay dlambda = (lambda - lambdal) / s;
800a7e14dcfSSatish Balay lambda_new = 0.75 * lambdau + 0.25 * lambda;
8019371c9d4SSatish Balay if (lambda_new > (lambda + dlambda)) lambda_new = lambda + dlambda;
802a7e14dcfSSatish Balay lambdal = lambda;
803a7e14dcfSSatish Balay rl = r;
804a7e14dcfSSatish Balay lambda = lambda_new;
805a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda);
806a7e14dcfSSatish Balay }
807a7e14dcfSSatish Balay }
808a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u);
809a7e14dcfSSatish Balay }
810a7e14dcfSSatish Balay
811a7e14dcfSSatish Balay *lam_ext = lambda;
8123387f462SBarry Smith if (innerIter >= df->maxProjIter) PetscCallAbort(PETSC_COMM_SELF, PetscInfo(NULL, "WARNING: DaiFletcher max iterations\n"));
813a7e14dcfSSatish Balay return innerIter;
814a7e14dcfSSatish Balay }
815