xref: /petsc/src/tao/leastsquares/impls/brgn/brgn.c (revision 30eeff36723c04aaacb4fcf78bbedf54c3f9917c)
1737f463aSAlp Dener #include <../src/tao/leastsquares/impls/brgn/brgn.h>
2737f463aSAlp Dener 
3a3c390cfSAlp Dener #define BRGN_user               0
4a3c390cfSAlp Dener #define BRGN_l2prox             1
5a3c390cfSAlp Dener #define BRGN_l1dict             2
6a3c390cfSAlp Dener #define BRGNRegTypes            3
7a3c390cfSAlp Dener 
8a3c390cfSAlp Dener static const char *BRGN_Table[64] = {"user","l2prox","l1dict"};
9a3c390cfSAlp Dener 
100d71dc2bSXiang Huang static PetscErrorCode GNHessianProd(Mat H,Vec in,Vec out)
110d71dc2bSXiang Huang {
120d71dc2bSXiang Huang   TAO_BRGN              *gn;
130d71dc2bSXiang Huang   PetscErrorCode        ierr;
140d71dc2bSXiang Huang 
150d71dc2bSXiang Huang   PetscFunctionBegin;
160d71dc2bSXiang Huang   ierr = MatShellGetContext(H,&gn);CHKERRQ(ierr);
170d71dc2bSXiang Huang   ierr = MatMult(gn->subsolver->ls_jac,in,gn->r_work);CHKERRQ(ierr);
180d71dc2bSXiang Huang   ierr = MatMultTranspose(gn->subsolver->ls_jac,gn->r_work,out);CHKERRQ(ierr);
19a3c390cfSAlp Dener   switch (gn->reg_type) {
20a3c390cfSAlp Dener   case BRGN_user:
21a3c390cfSAlp Dener     ierr = MatMult(gn->Hreg,in,gn->x_work);CHKERRQ(ierr);
228ac80d48SXiang Huang     ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr);
23a3c390cfSAlp Dener     break;
24a3c390cfSAlp Dener   case BRGN_l2prox:
25a3c390cfSAlp Dener     ierr = VecAXPY(out,gn->lambda,in);CHKERRQ(ierr);
26a3c390cfSAlp Dener     break;
27a3c390cfSAlp Dener   case BRGN_l1dict:
28a3c390cfSAlp Dener     /* out = out + lambda*D'*(diag.*(D*in)) */
29a3c390cfSAlp Dener     if (gn->D) {
30a3c390cfSAlp Dener       ierr = MatMult(gn->D,in,gn->y);CHKERRQ(ierr);/* y = D*in */
31a3c390cfSAlp Dener     } else {
32a3c390cfSAlp Dener       ierr = VecCopy(in,gn->y);CHKERRQ(ierr);
33a3c390cfSAlp Dener     }
34a3c390cfSAlp Dener     ierr = VecPointwiseMult(gn->y_work,gn->diag,gn->y);CHKERRQ(ierr);   /* y_work = diag.*(D*in), where diag = epsilon^2 ./ sqrt(x.^2+epsilon^2).^3 */
35a3c390cfSAlp Dener     if (gn->D) {
36a3c390cfSAlp Dener       ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr); /* x_work = D'*(diag.*(D*in)) */
37a3c390cfSAlp Dener     } else {
38a3c390cfSAlp Dener       ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr);
39a3c390cfSAlp Dener     }
40a3c390cfSAlp Dener     ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr);
41a3c390cfSAlp Dener     break;
42a3c390cfSAlp Dener   }
430d71dc2bSXiang Huang 
440d71dc2bSXiang Huang   PetscFunctionReturn(0);
450d71dc2bSXiang Huang }
460d71dc2bSXiang Huang 
470d71dc2bSXiang Huang static PetscErrorCode GNObjectiveGradientEval(Tao tao,Vec X,PetscReal *fcn,Vec G,void *ptr)
480d71dc2bSXiang Huang {
490d71dc2bSXiang Huang   TAO_BRGN              *gn = (TAO_BRGN *)ptr;
508e85b1b3SXiang Huang   PetscInt              K;                    /* dimension of D*X */
517cea06e1SXiang Huang   PetscScalar           yESum;
520d71dc2bSXiang Huang   PetscErrorCode        ierr;
53a3c390cfSAlp Dener   PetscReal             f_reg;
540d71dc2bSXiang Huang 
550d71dc2bSXiang Huang   PetscFunctionBegin;
568e85b1b3SXiang Huang     /* compute objective *fcn*/
57a3c390cfSAlp Dener   /* compute first term 0.5*||ls_res||_2^2 */
580d71dc2bSXiang Huang   ierr = TaoComputeResidual(tao,X,tao->ls_res);CHKERRQ(ierr);
59a3c390cfSAlp Dener   ierr = VecDot(tao->ls_res,tao->ls_res,fcn);CHKERRQ(ierr);
60a3c390cfSAlp Dener   *fcn *= 0.5;
61a3c390cfSAlp Dener   /* compute gradient G */
62a3c390cfSAlp Dener   ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr);
63a3c390cfSAlp Dener   ierr = MatMultTranspose(tao->ls_jac,tao->ls_res,G);CHKERRQ(ierr);
64a3c390cfSAlp Dener   /* add the regularization contribution */
65a3c390cfSAlp Dener   switch (gn->reg_type) {
66a3c390cfSAlp Dener   case BRGN_user:
67a3c390cfSAlp Dener     ierr = (*gn->regularizerobjandgrad)(tao,X,&f_reg,gn->x_work,gn->reg_obj_ctx);CHKERRQ(ierr);
68a3c390cfSAlp Dener     *fcn += gn->lambda*f_reg;
69a3c390cfSAlp Dener     ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr);
70a3c390cfSAlp Dener     break;
71a3c390cfSAlp Dener   case BRGN_l2prox:
72a3c390cfSAlp Dener     /* compute f = f + lambda*0.5*(xk - xkm1)^T(xk - xkm1) */
73a3c390cfSAlp Dener     ierr = VecAXPBYPCZ(gn->x_work,1.0,-1.0,0.0,X,gn->x_old);CHKERRQ(ierr); /*TODO: no need to use VecAXPBYPCZ for x - xkm1 */
74a3c390cfSAlp Dener     ierr = VecDot(gn->x_work,gn->x_work,&f_reg);CHKERRQ(ierr);
75a3c390cfSAlp Dener     *fcn += gn->lambda*0.5*f_reg;
76a3c390cfSAlp Dener     /* compute G = G + lambda*(xk - xkm1) */
77a3c390cfSAlp Dener     ierr = VecAXPBYPCZ(G,gn->lambda,-gn->lambda,1.0,X,gn->x_old);CHKERRQ(ierr);
78a3c390cfSAlp Dener     break;
79a3c390cfSAlp Dener   case BRGN_l1dict:
80a3c390cfSAlp Dener     /* compute f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x*/
81a3c390cfSAlp Dener     if (gn->D) {
827cea06e1SXiang Huang       ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */
83a3c390cfSAlp Dener     } else {
84a3c390cfSAlp Dener       ierr = VecCopy(X,gn->y);CHKERRQ(ierr);
85a3c390cfSAlp Dener     }
867cea06e1SXiang Huang     ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr);
877cea06e1SXiang Huang     ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
887cea06e1SXiang Huang     ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr);  /* gn->y_work = sqrt(y.^2+epsilon^2) */
897cea06e1SXiang Huang     ierr = VecSum(gn->y_work,&yESum);CHKERRQ(ierr);CHKERRQ(ierr);
908e85b1b3SXiang Huang     ierr = VecGetSize(gn->y,&K);CHKERRQ(ierr);
91a3c390cfSAlp Dener     *fcn += gn->lambda*(yESum - K*gn->epsilon);
927cea06e1SXiang Huang     /* compute G = G + lambda*D'*(y./sqrt(y.^2+epsilon^2)),where y = D*x */
937cea06e1SXiang Huang     ierr = VecPointwiseDivide(gn->y_work,gn->y,gn->y_work);CHKERRQ(ierr); /* reuse y_work = y./sqrt(y.^2+epsilon^2) */
94a3c390cfSAlp Dener     if (gn->D) {
957cea06e1SXiang Huang       ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr);
96a3c390cfSAlp Dener     } else {
97a3c390cfSAlp Dener       ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr);
98a3c390cfSAlp Dener     }
998ac80d48SXiang Huang     ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr);
100a3c390cfSAlp Dener     break;
101a3c390cfSAlp Dener   }
1020d71dc2bSXiang Huang   PetscFunctionReturn(0);
1030d71dc2bSXiang Huang }
1040d71dc2bSXiang Huang 
105737f463aSAlp Dener static PetscErrorCode GNComputeHessian(Tao tao,Vec X,Mat H,Mat Hpre,void *ptr)
106737f463aSAlp Dener {
1078ac80d48SXiang Huang   TAO_BRGN              *gn = (TAO_BRGN *)ptr;
108737f463aSAlp Dener   PetscErrorCode ierr;
109737f463aSAlp Dener 
110737f463aSAlp Dener   PetscFunctionBegin;
111e1e80dc8SAlp Dener   ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr);
1120d71dc2bSXiang Huang 
113a3c390cfSAlp Dener   switch (gn->reg_type) {
114a3c390cfSAlp Dener   case BRGN_user:
115a3c390cfSAlp Dener     ierr = (*gn->regularizerhessian)(tao,X,gn->Hreg,gn->reg_hess_ctx);CHKERRQ(ierr);
116a3c390cfSAlp Dener     break;
117a3c390cfSAlp Dener   case BRGN_l2prox:
118a3c390cfSAlp Dener     break;
119a3c390cfSAlp Dener   case BRGN_l1dict:
1207cea06e1SXiang Huang     /* calculate and store diagonal matrix as a vector: diag = epsilon^2 ./ sqrt(x.^2+epsilon^2).^3* --> diag = epsilon^2 ./ sqrt(y.^2+epsilon^2).^3,where y = D*x */
121a3c390cfSAlp Dener     if (gn->D) {
1227cea06e1SXiang Huang       ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */
123a3c390cfSAlp Dener     } else {
124a3c390cfSAlp Dener       ierr = VecCopy(X,gn->y);CHKERRQ(ierr);
125a3c390cfSAlp Dener     }
1267cea06e1SXiang Huang     ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr);
1277cea06e1SXiang Huang     ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
1287cea06e1SXiang Huang     ierr = VecCopy(gn->y_work,gn->diag);CHKERRQ(ierr);                  /* gn->diag = y.^2+epsilon^2 */
1297cea06e1SXiang Huang     ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr);                        /* gn->y_work = sqrt(y.^2+epsilon^2) */
1307cea06e1SXiang Huang     ierr = VecPointwiseMult(gn->diag,gn->y_work,gn->diag);CHKERRQ(ierr);/* gn->diag = sqrt(y.^2+epsilon^2).^3 */
1318ac80d48SXiang Huang     ierr = VecReciprocal(gn->diag);CHKERRQ(ierr);
1328ac80d48SXiang Huang     ierr = VecScale(gn->diag,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
133a3c390cfSAlp Dener     break;
134a3c390cfSAlp Dener   }
1358ac80d48SXiang Huang 
136e1e80dc8SAlp Dener   PetscFunctionReturn(0);
137e1e80dc8SAlp Dener }
138e1e80dc8SAlp Dener 
139e1e80dc8SAlp Dener static PetscErrorCode GNHookFunction(Tao tao,PetscInt iter)
140e1e80dc8SAlp Dener {
141e1e80dc8SAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->user_update;
142e1e80dc8SAlp Dener   PetscErrorCode        ierr;
143e1e80dc8SAlp Dener 
144e1e80dc8SAlp Dener   PetscFunctionBegin;
145e1e80dc8SAlp Dener   /* Update basic tao information from the subsolver */
146e1e80dc8SAlp Dener   gn->parent->nfuncs = tao->nfuncs;
147e1e80dc8SAlp Dener   gn->parent->ngrads = tao->ngrads;
148e1e80dc8SAlp Dener   gn->parent->nfuncgrads = tao->nfuncgrads;
149e1e80dc8SAlp Dener   gn->parent->nhess = tao->nhess;
150e1e80dc8SAlp Dener   gn->parent->niter = tao->niter;
151e1e80dc8SAlp Dener   gn->parent->ksp_its = tao->ksp_its;
152e1e80dc8SAlp Dener   gn->parent->ksp_tot_its = tao->ksp_tot_its;
153e1e80dc8SAlp Dener   ierr = TaoGetConvergedReason(tao,&gn->parent->reason);CHKERRQ(ierr);
154e1e80dc8SAlp Dener   /* Update the solution vectors */
155e1e80dc8SAlp Dener   if (iter == 0) {
156e1e80dc8SAlp Dener     ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr);
157e1e80dc8SAlp Dener   } else {
158e1e80dc8SAlp Dener     ierr = VecCopy(tao->solution,gn->x_old);CHKERRQ(ierr);
159e1e80dc8SAlp Dener     ierr = VecCopy(tao->solution,gn->parent->solution);CHKERRQ(ierr);
160e1e80dc8SAlp Dener   }
161e1e80dc8SAlp Dener   /* Update the gradient */
162e1e80dc8SAlp Dener   ierr = VecCopy(tao->gradient,gn->parent->gradient);CHKERRQ(ierr);
163e1e80dc8SAlp Dener   /* Call general purpose update function */
164e1e80dc8SAlp Dener   if (gn->parent->ops->update) {
165e1e80dc8SAlp Dener     ierr = (*gn->parent->ops->update)(gn->parent,gn->parent->niter);CHKERRQ(ierr);
166737f463aSAlp Dener   }
167737f463aSAlp Dener   PetscFunctionReturn(0);
168737f463aSAlp Dener }
169737f463aSAlp Dener 
170737f463aSAlp Dener static PetscErrorCode TaoSolve_BRGN(Tao tao)
171737f463aSAlp Dener {
172737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
173737f463aSAlp Dener   PetscErrorCode        ierr;
174737f463aSAlp Dener 
175737f463aSAlp Dener   PetscFunctionBegin;
176737f463aSAlp Dener   ierr = TaoSolve(gn->subsolver);CHKERRQ(ierr);
177e1e80dc8SAlp Dener   /* Update basic tao information from the subsolver */
178e1e80dc8SAlp Dener   tao->nfuncs = gn->subsolver->nfuncs;
179e1e80dc8SAlp Dener   tao->ngrads = gn->subsolver->ngrads;
180e1e80dc8SAlp Dener   tao->nfuncgrads = gn->subsolver->nfuncgrads;
181e1e80dc8SAlp Dener   tao->nhess = gn->subsolver->nhess;
182e1e80dc8SAlp Dener   tao->niter = gn->subsolver->niter;
183e1e80dc8SAlp Dener   tao->ksp_its = gn->subsolver->ksp_its;
184e1e80dc8SAlp Dener   tao->ksp_tot_its = gn->subsolver->ksp_tot_its;
185e1e80dc8SAlp Dener   ierr = TaoGetConvergedReason(gn->subsolver,&tao->reason);CHKERRQ(ierr);
186e1e80dc8SAlp Dener   /* Update vectors */
187e1e80dc8SAlp Dener   ierr = VecCopy(gn->subsolver->solution,tao->solution);CHKERRQ(ierr);
188e1e80dc8SAlp Dener   ierr = VecCopy(gn->subsolver->gradient,tao->gradient);CHKERRQ(ierr);
189737f463aSAlp Dener   PetscFunctionReturn(0);
190737f463aSAlp Dener }
191737f463aSAlp Dener 
192737f463aSAlp Dener static PetscErrorCode TaoSetFromOptions_BRGN(PetscOptionItems *PetscOptionsObject,Tao tao)
193737f463aSAlp Dener {
194737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
195737f463aSAlp Dener   PetscErrorCode        ierr;
196737f463aSAlp Dener 
197737f463aSAlp Dener   PetscFunctionBegin;
198a3c390cfSAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"least-squares problems with regularizer: ||f(x)||^2 + lambda*g(x), g(x) = ||xk-xkm1||^2 or ||Dx||_1 or user defined function.");CHKERRQ(ierr);
199*30eeff36SXiang Huang   ierr = PetscOptionsReal("-tao_brgn_lambda","regularizer weight (default 1e-4)","",gn->lambda,&gn->lambda,NULL);CHKERRQ(ierr);
200*30eeff36SXiang Huang   ierr = PetscOptionsReal("-tao_brgn_epsilon","L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6)","",gn->epsilon,&gn->epsilon,NULL);CHKERRQ(ierr);
201a3c390cfSAlp Dener   ierr = PetscOptionsEList("-tao_brgn_reg_type","regularization type", "",BRGN_Table,BRGNRegTypes,BRGN_Table[gn->reg_type],&gn->reg_type,NULL);CHKERRQ(ierr);
202737f463aSAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
203737f463aSAlp Dener   ierr = TaoSetFromOptions(gn->subsolver);CHKERRQ(ierr);
204737f463aSAlp Dener   PetscFunctionReturn(0);
205737f463aSAlp Dener }
206737f463aSAlp Dener 
207737f463aSAlp Dener static PetscErrorCode TaoView_BRGN(Tao tao,PetscViewer viewer)
208737f463aSAlp Dener {
209737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
210737f463aSAlp Dener   PetscErrorCode        ierr;
211737f463aSAlp Dener 
212737f463aSAlp Dener   PetscFunctionBegin;
213e1e80dc8SAlp Dener   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
214737f463aSAlp Dener   ierr = TaoView(gn->subsolver,viewer);CHKERRQ(ierr);
215e1e80dc8SAlp Dener   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
216737f463aSAlp Dener   PetscFunctionReturn(0);
217737f463aSAlp Dener }
218737f463aSAlp Dener 
219737f463aSAlp Dener static PetscErrorCode TaoSetUp_BRGN(Tao tao)
220737f463aSAlp Dener {
221737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
222737f463aSAlp Dener   PetscErrorCode        ierr;
223737f463aSAlp Dener   PetscBool             is_bnls,is_bntr,is_bntl;
2248e85b1b3SXiang Huang   PetscInt              i,n,N,K; /* dict has size K*N*/
225737f463aSAlp Dener 
226737f463aSAlp Dener   PetscFunctionBegin;
227737f463aSAlp Dener   if (!tao->ls_res) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualRoutine() must be called before setup!");
228737f463aSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNLS,&is_bnls);CHKERRQ(ierr);
229737f463aSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTR,&is_bntr);CHKERRQ(ierr);
230737f463aSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTL,&is_bntl);CHKERRQ(ierr);
231737f463aSAlp Dener   if ((is_bnls || is_bntr || is_bntl) && !tao->ls_jac) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualJacobianRoutine() must be called before setup!");
232e1e80dc8SAlp Dener   if (!tao->gradient) {
233e1e80dc8SAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
234e1e80dc8SAlp Dener   }
235737f463aSAlp Dener   if (!gn->x_work) {
236737f463aSAlp Dener     ierr = VecDuplicate(tao->solution,&gn->x_work);CHKERRQ(ierr);
237737f463aSAlp Dener   }
238737f463aSAlp Dener   if (!gn->r_work) {
239737f463aSAlp Dener     ierr = VecDuplicate(tao->ls_res,&gn->r_work);CHKERRQ(ierr);
240737f463aSAlp Dener   }
241e1e80dc8SAlp Dener   if (!gn->x_old) {
242e1e80dc8SAlp Dener     ierr = VecDuplicate(tao->solution,&gn->x_old);CHKERRQ(ierr);
243e1e80dc8SAlp Dener     ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr);
244e1e80dc8SAlp Dener   }
2457cea06e1SXiang Huang 
2468e85b1b3SXiang Huang   /*ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);*/
247*30eeff36SXiang Huang   if (BRGN_l1dict == gn->reg_type) {
248*30eeff36SXiang Huang     if (gn->D) {
249*30eeff36SXiang Huang       /* XH: debug: check matrix */
250*30eeff36SXiang Huang #if 0
251*30eeff36SXiang Huang       ierr = PetscPrintf(PETSC_COMM_SELF,"-------- Check D matrix: -------- \n"); CHKERRQ(ierr);
252*30eeff36SXiang Huang       ierr = MatView(gn->D,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
253*30eeff36SXiang Huang #endif
254*30eeff36SXiang Huang       ierr = MatGetSize(gn->D,&K,&N);CHKERRQ(ierr); /* Shell matrices still must have sizes defined. K = N for identity matrix, K=N-1 or N for gradient matrix */
255*30eeff36SXiang Huang     } else {
256*30eeff36SXiang Huang       ierr = VecGetSize(tao->solution,&K);CHKERRQ(ierr); /* If user does not setup dict matrix, use identiy matrix, K=N */
257*30eeff36SXiang Huang     }
2587cea06e1SXiang Huang     if (!gn->y) {
2597cea06e1SXiang Huang       ierr = VecCreate(PETSC_COMM_SELF,&gn->y);CHKERRQ(ierr);
2608e85b1b3SXiang Huang       ierr = VecSetSizes(gn->y,PETSC_DECIDE,K);CHKERRQ(ierr);
2617cea06e1SXiang Huang       ierr = VecSetFromOptions(gn->y);CHKERRQ(ierr);
2627cea06e1SXiang Huang       ierr = VecSet(gn->y,0.0);CHKERRQ(ierr);
2637cea06e1SXiang Huang 
2647cea06e1SXiang Huang     }
2657cea06e1SXiang Huang     if (!gn->y_work) {
2667cea06e1SXiang Huang       ierr = VecDuplicate(gn->y,&gn->y_work);CHKERRQ(ierr);
2677cea06e1SXiang Huang     }
2688ac80d48SXiang Huang     if (!gn->diag) {
2697cea06e1SXiang Huang       ierr = VecDuplicate(gn->y,&gn->diag);CHKERRQ(ierr);
2708ac80d48SXiang Huang       ierr = VecSet(gn->diag,0.0);CHKERRQ(ierr);
2718ac80d48SXiang Huang     }
272*30eeff36SXiang Huang   }
2730d71dc2bSXiang Huang 
2747cea06e1SXiang Huang 
275e1e80dc8SAlp Dener   if (!tao->setupcalled) {
276737f463aSAlp Dener     /* Hessian setup */
2778e85b1b3SXiang Huang     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
2788e85b1b3SXiang Huang     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
2798e85b1b3SXiang Huang     ierr = MatSetSizes(gn->H,n,n,N,N);CHKERRQ(ierr);
280737f463aSAlp Dener     ierr = MatSetType(gn->H,MATSHELL);CHKERRQ(ierr);
281737f463aSAlp Dener     ierr = MatSetUp(gn->H);CHKERRQ(ierr);
282737f463aSAlp Dener     ierr = MatShellSetOperation(gn->H,MATOP_MULT,(void (*)(void))GNHessianProd);CHKERRQ(ierr);
283737f463aSAlp Dener     ierr = MatShellSetContext(gn->H,(void*)gn);CHKERRQ(ierr);
2848e85b1b3SXiang Huang     /* Subsolver setup,include initial vector and dicttionary D */
285e1e80dc8SAlp Dener     ierr = TaoSetUpdate(gn->subsolver,GNHookFunction,(void*)gn);CHKERRQ(ierr);
286737f463aSAlp Dener     ierr = TaoSetInitialVector(gn->subsolver,tao->solution);CHKERRQ(ierr);
287737f463aSAlp Dener     if (tao->bounded) {
288737f463aSAlp Dener       ierr = TaoSetVariableBounds(gn->subsolver,tao->XL,tao->XU);CHKERRQ(ierr);
289737f463aSAlp Dener     }
290737f463aSAlp Dener     ierr = TaoSetResidualRoutine(gn->subsolver,tao->ls_res,tao->ops->computeresidual,tao->user_lsresP);CHKERRQ(ierr);
2914ffbe8acSAlp Dener     ierr = TaoSetJacobianResidualRoutine(gn->subsolver,tao->ls_jac,tao->ls_jac,tao->ops->computeresidualjacobian,tao->user_lsjacP);CHKERRQ(ierr);
292737f463aSAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(gn->subsolver,GNObjectiveGradientEval,(void*)gn);CHKERRQ(ierr);
293737f463aSAlp Dener     ierr = TaoSetHessianRoutine(gn->subsolver,gn->H,gn->H,GNComputeHessian,(void*)gn);CHKERRQ(ierr);
294e1e80dc8SAlp Dener     /* Propagate some options down */
295e1e80dc8SAlp Dener     ierr = TaoSetTolerances(gn->subsolver,tao->gatol,tao->grtol,tao->gttol);CHKERRQ(ierr);
296e1e80dc8SAlp Dener     ierr = TaoSetMaximumIterations(gn->subsolver,tao->max_it);CHKERRQ(ierr);
297e1e80dc8SAlp Dener     ierr = TaoSetMaximumFunctionEvaluations(gn->subsolver,tao->max_funcs);CHKERRQ(ierr);
298737f463aSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
299737f463aSAlp Dener       ierr = TaoSetMonitor(gn->subsolver,tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i]);CHKERRQ(ierr);
300737f463aSAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
301737f463aSAlp Dener     }
302737f463aSAlp Dener     ierr = TaoSetUp(gn->subsolver);CHKERRQ(ierr);
303e1e80dc8SAlp Dener   }
304737f463aSAlp Dener   PetscFunctionReturn(0);
305737f463aSAlp Dener }
306737f463aSAlp Dener 
307737f463aSAlp Dener static PetscErrorCode TaoDestroy_BRGN(Tao tao)
308737f463aSAlp Dener {
309737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
310737f463aSAlp Dener   PetscErrorCode        ierr;
311737f463aSAlp Dener 
312737f463aSAlp Dener   PetscFunctionBegin;
313737f463aSAlp Dener   if (tao->setupcalled) {
314e1e80dc8SAlp Dener     ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr);
315737f463aSAlp Dener     ierr = VecDestroy(&gn->x_work);CHKERRQ(ierr);
316737f463aSAlp Dener     ierr = VecDestroy(&gn->r_work);CHKERRQ(ierr);
317e1e80dc8SAlp Dener     ierr = VecDestroy(&gn->x_old);CHKERRQ(ierr);
3188ac80d48SXiang Huang     ierr = VecDestroy(&gn->diag);CHKERRQ(ierr);
3197cea06e1SXiang Huang     ierr = VecDestroy(&gn->y);CHKERRQ(ierr);
3207cea06e1SXiang Huang     ierr = VecDestroy(&gn->y_work);CHKERRQ(ierr);
321737f463aSAlp Dener   }
322737f463aSAlp Dener   ierr = MatDestroy(&gn->H);CHKERRQ(ierr);
3237cea06e1SXiang Huang   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
324737f463aSAlp Dener   ierr = TaoDestroy(&gn->subsolver);CHKERRQ(ierr);
325e1e80dc8SAlp Dener   gn->parent = NULL;
326737f463aSAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
327737f463aSAlp Dener   PetscFunctionReturn(0);
328737f463aSAlp Dener }
329737f463aSAlp Dener 
3303850be85SAlp Dener /*MC
3313850be85SAlp Dener   TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares
3323850be85SAlp Dener             problems with bound constraints. This algorithm is a thin wrapper around TAOBNTL
3333850be85SAlp Dener             that constructs the Guass-Newton problem with the user-provided least-squares
3343850be85SAlp Dener             residual and Jacobian. The problem is regularized with an L2-norm proximal point
3353850be85SAlp Dener             term.
3363850be85SAlp Dener 
3373850be85SAlp Dener   Options Database Keys:
3383850be85SAlp Dener   + -tao_bqnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
3393850be85SAlp Dener   . -tao_bqnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
3403850be85SAlp Dener   . -tao_bqnk_update_type - trust radius update method ("step", "direction", "interpolation")
3413850be85SAlp Dener   - -tao_bqnk_as_type - active-set estimation method ("none", "bertsekas")
3423850be85SAlp Dener 
3433850be85SAlp Dener   Level: beginner
3443850be85SAlp Dener M*/
345737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao)
346737f463aSAlp Dener {
347737f463aSAlp Dener   TAO_BRGN       *gn;
348737f463aSAlp Dener   PetscErrorCode ierr;
349737f463aSAlp Dener 
350737f463aSAlp Dener   PetscFunctionBegin;
351737f463aSAlp Dener   ierr = PetscNewLog(tao,&gn);CHKERRQ(ierr);
352737f463aSAlp Dener 
353737f463aSAlp Dener   tao->ops->destroy = TaoDestroy_BRGN;
354737f463aSAlp Dener   tao->ops->setup = TaoSetUp_BRGN;
355737f463aSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BRGN;
356737f463aSAlp Dener   tao->ops->view = TaoView_BRGN;
357737f463aSAlp Dener   tao->ops->solve = TaoSolve_BRGN;
358737f463aSAlp Dener 
359737f463aSAlp Dener   tao->data = (void*)gn;
360e1e80dc8SAlp Dener   gn->lambda = 1e-4;
3618ac80d48SXiang Huang   gn->epsilon = 1e-6;
362e1e80dc8SAlp Dener   gn->parent = tao;
363737f463aSAlp Dener 
364737f463aSAlp Dener   ierr = MatCreate(PetscObjectComm((PetscObject)tao),&gn->H);CHKERRQ(ierr);
365737f463aSAlp Dener   ierr = MatSetOptionsPrefix(gn->H,"tao_brgn_hessian_");CHKERRQ(ierr);
366737f463aSAlp Dener 
367737f463aSAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&gn->subsolver);CHKERRQ(ierr);
368737f463aSAlp Dener   ierr = TaoSetType(gn->subsolver,TAOBNLS);CHKERRQ(ierr);
369737f463aSAlp Dener   ierr = TaoSetOptionsPrefix(gn->subsolver,"tao_brgn_subsolver_");CHKERRQ(ierr);
370e1e80dc8SAlp Dener   PetscFunctionReturn(0);
371e1e80dc8SAlp Dener }
372e1e80dc8SAlp Dener 
373e1e80dc8SAlp Dener /*@C
374e1e80dc8SAlp Dener   TaoBRGNGetSubsolver - Get the pointer to the subsolver inside BRGN
375e1e80dc8SAlp Dener 
376e1e80dc8SAlp Dener   Collective on Tao
377e1e80dc8SAlp Dener 
378e1e80dc8SAlp Dener   Level: developer
379e1e80dc8SAlp Dener 
380e1e80dc8SAlp Dener   Input Parameters:
381e1e80dc8SAlp Dener +  tao - the Tao solver context
382e1e80dc8SAlp Dener -  subsolver - the Tao sub-solver context
383e1e80dc8SAlp Dener @*/
384e1e80dc8SAlp Dener PetscErrorCode TaoBRGNGetSubsolver(Tao tao,Tao *subsolver)
385e1e80dc8SAlp Dener {
386e1e80dc8SAlp Dener   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
387e1e80dc8SAlp Dener 
388e1e80dc8SAlp Dener   PetscFunctionBegin;
389e1e80dc8SAlp Dener   *subsolver = gn->subsolver;
390737f463aSAlp Dener   PetscFunctionReturn(0);
391737f463aSAlp Dener }
392737f463aSAlp Dener 
393737f463aSAlp Dener /*@C
3948e85b1b3SXiang Huang   TaoBRGNSetL1RegularizerWeight - Set the L1-norm regularizer weight for the Gauss-Newton least-squares algorithm
395737f463aSAlp Dener 
396737f463aSAlp Dener   Collective on Tao
397737f463aSAlp Dener 
398737f463aSAlp Dener   Level: developer
399737f463aSAlp Dener 
400737f463aSAlp Dener   Input Parameters:
401737f463aSAlp Dener +  tao - the Tao solver context
4028e85b1b3SXiang Huang -  lambda - L1-norm regularizer weight
403737f463aSAlp Dener @*/
404a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao,PetscReal lambda)
405737f463aSAlp Dener {
406737f463aSAlp Dener   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
407737f463aSAlp Dener 
4088ac80d48SXiang Huang   /* Initialize lambda here */
4090d71dc2bSXiang Huang 
410737f463aSAlp Dener   PetscFunctionBegin;
411737f463aSAlp Dener   gn->lambda = lambda;
412737f463aSAlp Dener   PetscFunctionReturn(0);
413737f463aSAlp Dener }
4140d71dc2bSXiang Huang 
4158ac80d48SXiang Huang /*@C
4168ac80d48SXiang Huang   TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm
4178ac80d48SXiang Huang 
4188ac80d48SXiang Huang   Collective on Tao
4198ac80d48SXiang Huang 
4208ac80d48SXiang Huang   Level: developer
4218ac80d48SXiang Huang 
4228ac80d48SXiang Huang   Input Parameters:
4238ac80d48SXiang Huang +  tao - the Tao solver context
4248ac80d48SXiang Huang -  epsilon - L1-norm smooth approximation parameter
4258ac80d48SXiang Huang @*/
4268ac80d48SXiang Huang PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao,PetscReal epsilon)
4278ac80d48SXiang Huang {
4288ac80d48SXiang Huang   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
4298ac80d48SXiang Huang 
4308ac80d48SXiang Huang   /* Initialize epsilon here */
4318ac80d48SXiang Huang 
4328ac80d48SXiang Huang   PetscFunctionBegin;
4338ac80d48SXiang Huang   gn->epsilon = epsilon;
4348ac80d48SXiang Huang   PetscFunctionReturn(0);
4358ac80d48SXiang Huang }
4368e85b1b3SXiang Huang 
4378e85b1b3SXiang Huang /*@C
4388e85b1b3SXiang Huang    TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem)
4398e85b1b3SXiang Huang 
4408e85b1b3SXiang Huang    Input Parameters:
4418e85b1b3SXiang Huang +  tao  - the Tao context
4428e85b1b3SXiang Huang .  dict - the user specified dictionary matrix
4438e85b1b3SXiang Huang 
4448e85b1b3SXiang Huang     Level: developer
4458e85b1b3SXiang Huang @*/
4468e85b1b3SXiang Huang PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao,Mat dict)
4478e85b1b3SXiang Huang {
4488e85b1b3SXiang Huang   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
4498e85b1b3SXiang Huang   PetscErrorCode ierr;
4508e85b1b3SXiang Huang   PetscFunctionBegin;
4518e85b1b3SXiang Huang   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
4528e85b1b3SXiang Huang   if (dict) {
4538e85b1b3SXiang Huang     PetscValidHeaderSpecific(dict,MAT_CLASSID,2);
454a3c390cfSAlp Dener     PetscCheckSameComm(tao,1,dict,2);
4558e85b1b3SXiang Huang     ierr = PetscObjectReference((PetscObject)dict);CHKERRQ(ierr);
4568e85b1b3SXiang Huang   }
4578e85b1b3SXiang Huang   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
4588e85b1b3SXiang Huang   gn->D = dict;  /* We allow to set a null dictionary, which means we just use default identity matrix? */
4598e85b1b3SXiang Huang   PetscFunctionReturn(0);
4608e85b1b3SXiang Huang }
4618e85b1b3SXiang Huang 
462a3c390cfSAlp Dener /*@C
463a3c390cfSAlp Dener @*/
464a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (*func)(Tao,Vec,PetscReal *,Vec,void*),void *ctx)
465a3c390cfSAlp Dener {
466a3c390cfSAlp Dener   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
467a3c390cfSAlp Dener 
468a3c390cfSAlp Dener   PetscFunctionBegin;
469a3c390cfSAlp Dener   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
470a3c390cfSAlp Dener   if (ctx) {
471a3c390cfSAlp Dener     gn->reg_obj_ctx = ctx;
472a3c390cfSAlp Dener   }
473a3c390cfSAlp Dener   if (func) {
474a3c390cfSAlp Dener     gn->regularizerobjandgrad = func;
475a3c390cfSAlp Dener   }
476a3c390cfSAlp Dener   PetscFunctionReturn(0);
477a3c390cfSAlp Dener }
478a3c390cfSAlp Dener 
479a3c390cfSAlp Dener /*@C
480a3c390cfSAlp Dener @*/
481a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao,Mat Hreg,PetscErrorCode (*func)(Tao,Vec,Mat,void*),void *ctx)
482a3c390cfSAlp Dener {
483a3c390cfSAlp Dener   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
484a3c390cfSAlp Dener   PetscErrorCode ierr;
485a3c390cfSAlp Dener 
486a3c390cfSAlp Dener   PetscFunctionBegin;
487a3c390cfSAlp Dener   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
488a3c390cfSAlp Dener   if (Hreg) {
489a3c390cfSAlp Dener     PetscValidHeaderSpecific(Hreg,MAT_CLASSID,2);
490a3c390cfSAlp Dener     PetscCheckSameComm(tao,1,Hreg,2);
491a3c390cfSAlp Dener   } else {
492a3c390cfSAlp Dener     SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONG,"NULL Hessian detected! User must provide valid Hessian for the regularizer.");
493a3c390cfSAlp Dener   }
494a3c390cfSAlp Dener   if (ctx) {
495a3c390cfSAlp Dener     gn->reg_hess_ctx = ctx;
496a3c390cfSAlp Dener   }
497a3c390cfSAlp Dener   if (func) {
498a3c390cfSAlp Dener     gn->regularizerhessian = func;
499a3c390cfSAlp Dener   }
500a3c390cfSAlp Dener   if (Hreg) {
501a3c390cfSAlp Dener     ierr = PetscObjectReference((PetscObject)Hreg);CHKERRQ(ierr);
502a3c390cfSAlp Dener     ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr);
503a3c390cfSAlp Dener     gn->Hreg = Hreg;
504a3c390cfSAlp Dener   }
505a3c390cfSAlp Dener   PetscFunctionReturn(0);
506a3c390cfSAlp Dener }