xref: /petsc/src/tao/leastsquares/impls/brgn/brgn.c (revision 8e85b1b3a48c719bb7d116aee5cfc4a5b180d411)
1737f463aSAlp Dener #include <../src/tao/leastsquares/impls/brgn/brgn.h>
2737f463aSAlp Dener 
30d71dc2bSXiang Huang static PetscErrorCode GNHessianProd(Mat H, Vec in, Vec out)
40d71dc2bSXiang Huang {
50d71dc2bSXiang Huang   TAO_BRGN              *gn;
60d71dc2bSXiang Huang   PetscErrorCode        ierr;
70d71dc2bSXiang Huang 
80d71dc2bSXiang Huang   PetscFunctionBegin;
90d71dc2bSXiang Huang   ierr = MatShellGetContext(H, &gn);CHKERRQ(ierr);
100d71dc2bSXiang Huang   ierr = MatMult(gn->subsolver->ls_jac, in, gn->r_work);CHKERRQ(ierr);
110d71dc2bSXiang Huang   ierr = MatMultTranspose(gn->subsolver->ls_jac, gn->r_work, out);CHKERRQ(ierr);
127cea06e1SXiang Huang   /* out = out + lambda*D'*(diag.*(D*in)) */
137cea06e1SXiang Huang   ierr = MatMult(gn->D, in, gn->y);CHKERRQ(ierr);  /* y = D*in */
147cea06e1SXiang Huang   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 */
157cea06e1SXiang Huang   ierr = MatMultTranspose(gn->D, gn->y_work, gn->x_work);CHKERRQ(ierr);   /* x_work = D'*(diag.*(D*in)) */
168ac80d48SXiang Huang   ierr = VecAXPY(out, gn->lambda, gn->x_work);CHKERRQ(ierr);
170d71dc2bSXiang Huang 
180d71dc2bSXiang Huang   PetscFunctionReturn(0);
190d71dc2bSXiang Huang }
200d71dc2bSXiang Huang 
210d71dc2bSXiang Huang static PetscErrorCode GNObjectiveGradientEval(Tao tao, Vec X, PetscReal *fcn, Vec G, void *ptr)
220d71dc2bSXiang Huang {
230d71dc2bSXiang Huang   TAO_BRGN              *gn = (TAO_BRGN *)ptr;
24*8e85b1b3SXiang Huang   PetscInt              K;                    /* dimension of D*X */
257cea06e1SXiang Huang   PetscScalar           yESum;
260d71dc2bSXiang Huang   PetscErrorCode        ierr;
270d71dc2bSXiang Huang 
280d71dc2bSXiang Huang   PetscFunctionBegin;
29*8e85b1b3SXiang Huang     /* compute objective *fcn*/
308ac80d48SXiang Huang   /* compute first term ||ls_res||^2 */
310d71dc2bSXiang Huang   ierr = TaoComputeResidual(tao, X, tao->ls_res);CHKERRQ(ierr);
320d71dc2bSXiang Huang   ierr = VecDotBegin(tao->ls_res, tao->ls_res, fcn);CHKERRQ(ierr);
330d71dc2bSXiang Huang   ierr = VecDotEnd(tao->ls_res, tao->ls_res, fcn);CHKERRQ(ierr);
347cea06e1SXiang Huang   /* add the second term lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x*/
357cea06e1SXiang Huang   ierr = MatMult(gn->D, X, gn->y);CHKERRQ(ierr);  /* y = D*x */
367cea06e1SXiang Huang   ierr = VecPointwiseMult(gn->y_work, gn->y, gn->y);CHKERRQ(ierr);
377cea06e1SXiang Huang   ierr = VecShift(gn->y_work, gn->epsilon*gn->epsilon);CHKERRQ(ierr);
387cea06e1SXiang Huang   ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr);      /* gn->y_work = sqrt(y.^2+epsilon^2) */
397cea06e1SXiang Huang   ierr = VecSum(gn->y_work, &yESum);CHKERRQ(ierr);CHKERRQ(ierr);
40*8e85b1b3SXiang Huang   ierr = VecGetSize(gn->y, &K);CHKERRQ(ierr);
41*8e85b1b3SXiang Huang   *fcn = 0.5*(*fcn) + gn->lambda*(yESum - K*gn->epsilon);
420d71dc2bSXiang Huang 
438ac80d48SXiang Huang   /* compute gradient G */
440d71dc2bSXiang Huang   ierr = TaoComputeResidualJacobian(tao, X, tao->ls_jac, tao->ls_jac_pre);CHKERRQ(ierr);
450d71dc2bSXiang Huang   ierr = MatMultTranspose(tao->ls_jac, tao->ls_res, G);CHKERRQ(ierr);
467cea06e1SXiang Huang   /* compute G = G + lambda*D'*(y./sqrt(y.^2+epsilon^2)), where y = D*x */
477cea06e1SXiang Huang   ierr = VecPointwiseDivide(gn->y_work, gn->y, gn->y_work);CHKERRQ(ierr); /* reuse y_work = y./sqrt(y.^2+epsilon^2) */
487cea06e1SXiang Huang   ierr = MatMultTranspose(gn->D, gn->y_work, gn->x_work);CHKERRQ(ierr);
498ac80d48SXiang Huang   ierr = VecAXPY(G, gn->lambda, gn->x_work);CHKERRQ(ierr);
500d71dc2bSXiang Huang 
510d71dc2bSXiang Huang   PetscFunctionReturn(0);
520d71dc2bSXiang Huang }
530d71dc2bSXiang Huang 
54737f463aSAlp Dener 
55737f463aSAlp Dener static PetscErrorCode GNComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
56737f463aSAlp Dener {
578ac80d48SXiang Huang   TAO_BRGN              *gn = (TAO_BRGN *)ptr;
58737f463aSAlp Dener   PetscErrorCode ierr;
59737f463aSAlp Dener 
60737f463aSAlp Dener   PetscFunctionBegin;
61e1e80dc8SAlp Dener   ierr = TaoComputeResidualJacobian(tao, X, tao->ls_jac, tao->ls_jac_pre);CHKERRQ(ierr);
620d71dc2bSXiang Huang 
637cea06e1SXiang 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 */
647cea06e1SXiang Huang   ierr = MatMult(gn->D, X, gn->y);CHKERRQ(ierr);  /* y = D*x */
657cea06e1SXiang Huang   ierr = VecPointwiseMult(gn->y_work, gn->y, gn->y);CHKERRQ(ierr);
667cea06e1SXiang Huang   ierr = VecShift(gn->y_work, gn->epsilon*gn->epsilon);CHKERRQ(ierr);
677cea06e1SXiang Huang   ierr = VecCopy(gn->y_work, gn->diag);CHKERRQ(ierr);                     /* gn->diag = y.^2+epsilon^2 */
687cea06e1SXiang Huang   ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr);                            /* gn->y_work = sqrt(y.^2+epsilon^2) */
697cea06e1SXiang Huang   ierr = VecPointwiseMult(gn->diag, gn->y_work, gn->diag);CHKERRQ(ierr);  /* gn->diag = sqrt(y.^2+epsilon^2).^3 */
708ac80d48SXiang Huang   ierr = VecReciprocal(gn->diag);CHKERRQ(ierr);
718ac80d48SXiang Huang   ierr = VecScale(gn->diag, gn->epsilon*gn->epsilon);CHKERRQ(ierr);
728ac80d48SXiang Huang 
73e1e80dc8SAlp Dener   PetscFunctionReturn(0);
74e1e80dc8SAlp Dener }
75e1e80dc8SAlp Dener 
76e1e80dc8SAlp Dener static PetscErrorCode GNHookFunction(Tao tao, PetscInt iter)
77e1e80dc8SAlp Dener {
78e1e80dc8SAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->user_update;
79e1e80dc8SAlp Dener   PetscErrorCode        ierr;
80e1e80dc8SAlp Dener 
81e1e80dc8SAlp Dener   PetscFunctionBegin;
82e1e80dc8SAlp Dener   /* Update basic tao information from the subsolver */
83e1e80dc8SAlp Dener   gn->parent->nfuncs = tao->nfuncs;
84e1e80dc8SAlp Dener   gn->parent->ngrads = tao->ngrads;
85e1e80dc8SAlp Dener   gn->parent->nfuncgrads = tao->nfuncgrads;
86e1e80dc8SAlp Dener   gn->parent->nhess = tao->nhess;
87e1e80dc8SAlp Dener   gn->parent->niter = tao->niter;
88e1e80dc8SAlp Dener   gn->parent->ksp_its = tao->ksp_its;
89e1e80dc8SAlp Dener   gn->parent->ksp_tot_its = tao->ksp_tot_its;
90e1e80dc8SAlp Dener   ierr = TaoGetConvergedReason(tao, &gn->parent->reason);CHKERRQ(ierr);
91e1e80dc8SAlp Dener   /* Update the solution vectors */
92e1e80dc8SAlp Dener   if (iter == 0) {
93e1e80dc8SAlp Dener     ierr = VecSet(gn->x_old, 0.0);CHKERRQ(ierr);
94e1e80dc8SAlp Dener   } else {
95e1e80dc8SAlp Dener     ierr = VecCopy(tao->solution, gn->x_old);CHKERRQ(ierr);
96e1e80dc8SAlp Dener     ierr = VecCopy(tao->solution, gn->parent->solution);CHKERRQ(ierr);
97e1e80dc8SAlp Dener   }
98e1e80dc8SAlp Dener   /* Update the gradient */
99e1e80dc8SAlp Dener   ierr = VecCopy(tao->gradient, gn->parent->gradient);CHKERRQ(ierr);
100e1e80dc8SAlp Dener   /* Call general purpose update function */
101e1e80dc8SAlp Dener   if (gn->parent->ops->update) {
102e1e80dc8SAlp Dener     ierr = (*gn->parent->ops->update)(gn->parent, gn->parent->niter);CHKERRQ(ierr);
103737f463aSAlp Dener   }
104737f463aSAlp Dener   PetscFunctionReturn(0);
105737f463aSAlp Dener }
106737f463aSAlp Dener 
107737f463aSAlp Dener static PetscErrorCode TaoSolve_BRGN(Tao tao)
108737f463aSAlp Dener {
109737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
110737f463aSAlp Dener   PetscErrorCode        ierr;
111737f463aSAlp Dener 
112737f463aSAlp Dener   PetscFunctionBegin;
113737f463aSAlp Dener   ierr = TaoSolve(gn->subsolver);CHKERRQ(ierr);
114e1e80dc8SAlp Dener   /* Update basic tao information from the subsolver */
115e1e80dc8SAlp Dener   tao->nfuncs = gn->subsolver->nfuncs;
116e1e80dc8SAlp Dener   tao->ngrads = gn->subsolver->ngrads;
117e1e80dc8SAlp Dener   tao->nfuncgrads = gn->subsolver->nfuncgrads;
118e1e80dc8SAlp Dener   tao->nhess = gn->subsolver->nhess;
119e1e80dc8SAlp Dener   tao->niter = gn->subsolver->niter;
120e1e80dc8SAlp Dener   tao->ksp_its = gn->subsolver->ksp_its;
121e1e80dc8SAlp Dener   tao->ksp_tot_its = gn->subsolver->ksp_tot_its;
122e1e80dc8SAlp Dener   ierr = TaoGetConvergedReason(gn->subsolver, &tao->reason);CHKERRQ(ierr);
123e1e80dc8SAlp Dener   /* Update vectors */
124e1e80dc8SAlp Dener   ierr = VecCopy(gn->subsolver->solution, tao->solution);CHKERRQ(ierr);
125e1e80dc8SAlp Dener   ierr = VecCopy(gn->subsolver->gradient, tao->gradient);CHKERRQ(ierr);
126737f463aSAlp Dener   PetscFunctionReturn(0);
127737f463aSAlp Dener }
128737f463aSAlp Dener 
129737f463aSAlp Dener static PetscErrorCode TaoSetFromOptions_BRGN(PetscOptionItems *PetscOptionsObject,Tao tao)
130737f463aSAlp Dener {
131737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
132737f463aSAlp Dener   PetscErrorCode        ierr;
133737f463aSAlp Dener 
134737f463aSAlp Dener   PetscFunctionBegin;
1358ac80d48SXiang Huang   /* old Tikhonov regularization code
136737f463aSAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Gauss-Newton method for least-squares problems using Tikhonov regularization");CHKERRQ(ierr);
137737f463aSAlp Dener   ierr = PetscOptionsReal("-tao_brgn_lambda", "Tikhonov regularization factor", "", gn->lambda, &gn->lambda, NULL);CHKERRQ(ierr);
1388ac80d48SXiang Huang   */
1398ac80d48SXiang Huang   ierr = PetscOptionsHead(PetscOptionsObject,"least-squares problems with L1 regularizer: ||f(x)||^2 + lambda*||x||_1. Currently L1-norm is approximated with smooth form");CHKERRQ(ierr);
1408ac80d48SXiang Huang   ierr = PetscOptionsReal("-tao_brgn_lambda", "L1-norm regularizer weight", "", gn->lambda, &gn->lambda, NULL);CHKERRQ(ierr);
1418ac80d48SXiang Huang   ierr = PetscOptionsReal("-tao_brgn_epsilon", "L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon)", "", gn->epsilon, &gn->epsilon, NULL);CHKERRQ(ierr);
142737f463aSAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
143737f463aSAlp Dener   ierr = TaoSetFromOptions(gn->subsolver);CHKERRQ(ierr);
144737f463aSAlp Dener   PetscFunctionReturn(0);
145737f463aSAlp Dener }
146737f463aSAlp Dener 
147737f463aSAlp Dener static PetscErrorCode TaoView_BRGN(Tao tao, PetscViewer viewer)
148737f463aSAlp Dener {
149737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
150737f463aSAlp Dener   PetscErrorCode        ierr;
151737f463aSAlp Dener 
152737f463aSAlp Dener   PetscFunctionBegin;
153e1e80dc8SAlp Dener   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
154737f463aSAlp Dener   ierr = TaoView(gn->subsolver, viewer);CHKERRQ(ierr);
155e1e80dc8SAlp Dener   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
156737f463aSAlp Dener   PetscFunctionReturn(0);
157737f463aSAlp Dener }
158737f463aSAlp Dener 
159737f463aSAlp Dener static PetscErrorCode TaoSetUp_BRGN(Tao tao)
160737f463aSAlp Dener {
161737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
162737f463aSAlp Dener   PetscErrorCode        ierr;
163737f463aSAlp Dener   PetscBool             is_bnls, is_bntr, is_bntl;
164*8e85b1b3SXiang Huang   PetscInt              i, n, N, K; /* dict has size K*N*/
165*8e85b1b3SXiang Huang   /*PetscScalar           v; */ /* XH: hack to set value of matrix */
166737f463aSAlp Dener 
167737f463aSAlp Dener   PetscFunctionBegin;
168737f463aSAlp Dener   if (!tao->ls_res) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetResidualRoutine() must be called before setup!");
169737f463aSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNLS, &is_bnls);CHKERRQ(ierr);
170737f463aSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNTR, &is_bntr);CHKERRQ(ierr);
171737f463aSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNTL, &is_bntl);CHKERRQ(ierr);
172737f463aSAlp Dener   if ((is_bnls || is_bntr || is_bntl) && !tao->ls_jac) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetResidualJacobianRoutine() must be called before setup!");
173e1e80dc8SAlp Dener   if (!tao->gradient){
174e1e80dc8SAlp Dener     ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);
175e1e80dc8SAlp Dener   }
176737f463aSAlp Dener   if (!gn->x_work){
177737f463aSAlp Dener     ierr = VecDuplicate(tao->solution, &gn->x_work);CHKERRQ(ierr);
178737f463aSAlp Dener   }
179737f463aSAlp Dener   if (!gn->r_work){
180737f463aSAlp Dener     ierr = VecDuplicate(tao->ls_res, &gn->r_work);CHKERRQ(ierr);
181737f463aSAlp Dener   }
182e1e80dc8SAlp Dener   if (!gn->x_old) {
183e1e80dc8SAlp Dener     ierr = VecDuplicate(tao->solution, &gn->x_old);CHKERRQ(ierr);
184e1e80dc8SAlp Dener     ierr = VecSet(gn->x_old, 0.0);CHKERRQ(ierr);
185e1e80dc8SAlp Dener   }
1867cea06e1SXiang Huang 
187*8e85b1b3SXiang Huang   /*ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);*/
188*8e85b1b3SXiang Huang   /* TODO: Safeguard against NULL matrix */
189*8e85b1b3SXiang Huang   /*if (!gn->D)*/
190*8e85b1b3SXiang Huang   ierr = MatGetSize(gn->D, &K, &N);CHKERRQ(ierr); /* Shell matrices still must have sizes defined */
191*8e85b1b3SXiang Huang    /* K = N for identity matrix, K=N-1 or N for gradient matrix */
1927cea06e1SXiang Huang   if (!gn->y){
1937cea06e1SXiang Huang     ierr = VecCreate(PETSC_COMM_SELF,&gn->y);CHKERRQ(ierr);
194*8e85b1b3SXiang Huang     ierr = VecSetSizes(gn->y,PETSC_DECIDE,K);CHKERRQ(ierr);
1957cea06e1SXiang Huang     ierr = VecSetFromOptions(gn->y);CHKERRQ(ierr);
1967cea06e1SXiang Huang     ierr = VecSet(gn->y,0.0);CHKERRQ(ierr);
1977cea06e1SXiang Huang 
1987cea06e1SXiang Huang   }
1997cea06e1SXiang Huang   if (!gn->y_work){
2007cea06e1SXiang Huang     ierr = VecDuplicate(gn->y,&gn->y_work);CHKERRQ(ierr);
2017cea06e1SXiang Huang   }
2028ac80d48SXiang Huang   if (!gn->diag){
2037cea06e1SXiang Huang     ierr = VecDuplicate(gn->y,&gn->diag);CHKERRQ(ierr);
2048ac80d48SXiang Huang     ierr = VecSet(gn->diag,0.0);CHKERRQ(ierr);
2058ac80d48SXiang Huang   }
2060d71dc2bSXiang Huang 
2077cea06e1SXiang Huang   /* XH: debug: check matrix */
208*8e85b1b3SXiang Huang   ierr = PetscPrintf(PETSC_COMM_SELF, "-------- Check D matrix: -------- \n"); CHKERRQ(ierr);
2097cea06e1SXiang Huang   ierr = MatView(gn->D,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
21020fe612cSXiang Huang 
2117cea06e1SXiang Huang 
212e1e80dc8SAlp Dener   if (!tao->setupcalled) {
213737f463aSAlp Dener     /* Hessian setup */
214*8e85b1b3SXiang Huang     ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr);
215*8e85b1b3SXiang Huang     ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);
216*8e85b1b3SXiang Huang     ierr = MatSetSizes(gn->H, n, n, N, N);CHKERRQ(ierr);
217737f463aSAlp Dener     ierr = MatSetType(gn->H, MATSHELL);CHKERRQ(ierr);
218737f463aSAlp Dener     ierr = MatSetUp(gn->H);CHKERRQ(ierr);
219737f463aSAlp Dener     ierr = MatShellSetOperation(gn->H, MATOP_MULT, (void (*)(void))GNHessianProd);CHKERRQ(ierr);
220737f463aSAlp Dener     ierr = MatShellSetContext(gn->H, (void*)gn);CHKERRQ(ierr);
221*8e85b1b3SXiang Huang     /* Subsolver setup, include initial vector and dicttionary D */
222e1e80dc8SAlp Dener     ierr = TaoSetUpdate(gn->subsolver, GNHookFunction, (void*)gn);CHKERRQ(ierr);
223737f463aSAlp Dener     ierr = TaoSetInitialVector(gn->subsolver, tao->solution);CHKERRQ(ierr);
224737f463aSAlp Dener     if (tao->bounded) {
225737f463aSAlp Dener       ierr = TaoSetVariableBounds(gn->subsolver, tao->XL, tao->XU);CHKERRQ(ierr);
226737f463aSAlp Dener     }
227737f463aSAlp Dener     ierr = TaoSetResidualRoutine(gn->subsolver, tao->ls_res, tao->ops->computeresidual, tao->user_lsresP);CHKERRQ(ierr);
2284ffbe8acSAlp Dener     ierr = TaoSetJacobianResidualRoutine(gn->subsolver, tao->ls_jac, tao->ls_jac, tao->ops->computeresidualjacobian, tao->user_lsjacP);CHKERRQ(ierr);
229737f463aSAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(gn->subsolver, GNObjectiveGradientEval, (void*)gn);CHKERRQ(ierr);
230737f463aSAlp Dener     ierr = TaoSetHessianRoutine(gn->subsolver, gn->H, gn->H, GNComputeHessian, (void*)gn);CHKERRQ(ierr);
231e1e80dc8SAlp Dener     /* Propagate some options down */
232e1e80dc8SAlp Dener     ierr = TaoSetTolerances(gn->subsolver, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
233e1e80dc8SAlp Dener     ierr = TaoSetMaximumIterations(gn->subsolver, tao->max_it);CHKERRQ(ierr);
234e1e80dc8SAlp Dener     ierr = TaoSetMaximumFunctionEvaluations(gn->subsolver, tao->max_funcs);CHKERRQ(ierr);
235737f463aSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
236737f463aSAlp Dener       ierr = TaoSetMonitor(gn->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
237737f463aSAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
238737f463aSAlp Dener     }
239737f463aSAlp Dener     ierr = TaoSetUp(gn->subsolver);CHKERRQ(ierr);
240e1e80dc8SAlp Dener   }
241737f463aSAlp Dener   PetscFunctionReturn(0);
242737f463aSAlp Dener }
243737f463aSAlp Dener 
244737f463aSAlp Dener static PetscErrorCode TaoDestroy_BRGN(Tao tao)
245737f463aSAlp Dener {
246737f463aSAlp Dener   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
247737f463aSAlp Dener   PetscErrorCode        ierr;
248737f463aSAlp Dener 
249737f463aSAlp Dener   PetscFunctionBegin;
250737f463aSAlp Dener   if (tao->setupcalled) {
251e1e80dc8SAlp Dener     ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr);
252737f463aSAlp Dener     ierr = VecDestroy(&gn->x_work);CHKERRQ(ierr);
253737f463aSAlp Dener     ierr = VecDestroy(&gn->r_work);CHKERRQ(ierr);
254e1e80dc8SAlp Dener     ierr = VecDestroy(&gn->x_old);CHKERRQ(ierr);
2558ac80d48SXiang Huang     ierr = VecDestroy(&gn->diag);CHKERRQ(ierr);
2567cea06e1SXiang Huang     ierr = VecDestroy(&gn->y);CHKERRQ(ierr);
2577cea06e1SXiang Huang     ierr = VecDestroy(&gn->y_work);CHKERRQ(ierr);
258737f463aSAlp Dener   }
259737f463aSAlp Dener   ierr = MatDestroy(&gn->H);CHKERRQ(ierr);
2607cea06e1SXiang Huang   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
261737f463aSAlp Dener   ierr = TaoDestroy(&gn->subsolver);CHKERRQ(ierr);
262e1e80dc8SAlp Dener   gn->parent = NULL;
263737f463aSAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
264737f463aSAlp Dener   PetscFunctionReturn(0);
265737f463aSAlp Dener }
266737f463aSAlp Dener 
2673850be85SAlp Dener /*MC
2683850be85SAlp Dener   TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares
2693850be85SAlp Dener             problems with bound constraints. This algorithm is a thin wrapper around TAOBNTL
2703850be85SAlp Dener             that constructs the Guass-Newton problem with the user-provided least-squares
2713850be85SAlp Dener             residual and Jacobian. The problem is regularized with an L2-norm proximal point
2723850be85SAlp Dener             term.
2733850be85SAlp Dener 
2743850be85SAlp Dener   Options Database Keys:
2753850be85SAlp Dener   + -tao_bqnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
2763850be85SAlp Dener   . -tao_bqnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
2773850be85SAlp Dener   . -tao_bqnk_update_type - trust radius update method ("step", "direction", "interpolation")
2783850be85SAlp Dener   - -tao_bqnk_as_type - active-set estimation method ("none", "bertsekas")
2793850be85SAlp Dener 
2803850be85SAlp Dener   Level: beginner
2813850be85SAlp Dener M*/
282737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao)
283737f463aSAlp Dener {
284737f463aSAlp Dener   TAO_BRGN       *gn;
285737f463aSAlp Dener   PetscErrorCode ierr;
286737f463aSAlp Dener 
287737f463aSAlp Dener   PetscFunctionBegin;
288737f463aSAlp Dener   ierr = PetscNewLog(tao,&gn);CHKERRQ(ierr);
289737f463aSAlp Dener 
290737f463aSAlp Dener   tao->ops->destroy = TaoDestroy_BRGN;
291737f463aSAlp Dener   tao->ops->setup = TaoSetUp_BRGN;
292737f463aSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BRGN;
293737f463aSAlp Dener   tao->ops->view = TaoView_BRGN;
294737f463aSAlp Dener   tao->ops->solve = TaoSolve_BRGN;
295737f463aSAlp Dener 
296737f463aSAlp Dener   tao->data = (void*)gn;
297e1e80dc8SAlp Dener   gn->lambda = 1e-4;
2988ac80d48SXiang Huang   gn->epsilon = 1e-6;
299e1e80dc8SAlp Dener   gn->parent = tao;
300737f463aSAlp Dener 
301737f463aSAlp Dener   ierr = MatCreate(PetscObjectComm((PetscObject)tao), &gn->H);CHKERRQ(ierr);
302737f463aSAlp Dener   ierr = MatSetOptionsPrefix(gn->H, "tao_brgn_hessian_");CHKERRQ(ierr);
303737f463aSAlp Dener 
304737f463aSAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &gn->subsolver);CHKERRQ(ierr);
305737f463aSAlp Dener   ierr = TaoSetType(gn->subsolver, TAOBNLS);CHKERRQ(ierr);
306737f463aSAlp Dener   ierr = TaoSetOptionsPrefix(gn->subsolver, "tao_brgn_subsolver_");CHKERRQ(ierr);
307e1e80dc8SAlp Dener   PetscFunctionReturn(0);
308e1e80dc8SAlp Dener }
309e1e80dc8SAlp Dener 
310e1e80dc8SAlp Dener /*@C
311e1e80dc8SAlp Dener   TaoBRGNGetSubsolver - Get the pointer to the subsolver inside BRGN
312e1e80dc8SAlp Dener 
313e1e80dc8SAlp Dener   Collective on Tao
314e1e80dc8SAlp Dener 
315e1e80dc8SAlp Dener   Level: developer
316e1e80dc8SAlp Dener 
317e1e80dc8SAlp Dener   Input Parameters:
318e1e80dc8SAlp Dener +  tao - the Tao solver context
319e1e80dc8SAlp Dener -  subsolver - the Tao sub-solver context
320e1e80dc8SAlp Dener @*/
321e1e80dc8SAlp Dener PetscErrorCode TaoBRGNGetSubsolver(Tao tao,Tao *subsolver)
322e1e80dc8SAlp Dener {
323e1e80dc8SAlp Dener   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
324e1e80dc8SAlp Dener 
325e1e80dc8SAlp Dener   PetscFunctionBegin;
326e1e80dc8SAlp Dener   *subsolver = gn->subsolver;
327737f463aSAlp Dener   PetscFunctionReturn(0);
328737f463aSAlp Dener }
329737f463aSAlp Dener 
330737f463aSAlp Dener /*@C
331*8e85b1b3SXiang Huang   TaoBRGNSetL1RegularizerWeight - Set the L1-norm regularizer weight for the Gauss-Newton least-squares algorithm
332737f463aSAlp Dener 
333737f463aSAlp Dener   Collective on Tao
334737f463aSAlp Dener 
335737f463aSAlp Dener   Level: developer
336737f463aSAlp Dener 
337737f463aSAlp Dener   Input Parameters:
338737f463aSAlp Dener +  tao - the Tao solver context
339*8e85b1b3SXiang Huang -  lambda - L1-norm regularizer weight
340737f463aSAlp Dener @*/
341*8e85b1b3SXiang Huang PetscErrorCode TaoBRGNSetL1RegularizerWeight(Tao tao,PetscReal lambda)
342737f463aSAlp Dener {
343737f463aSAlp Dener   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
344737f463aSAlp Dener 
3458ac80d48SXiang Huang   /* Initialize lambda here */
3460d71dc2bSXiang Huang 
347737f463aSAlp Dener   PetscFunctionBegin;
348737f463aSAlp Dener   gn->lambda = lambda;
349737f463aSAlp Dener   PetscFunctionReturn(0);
350737f463aSAlp Dener }
3510d71dc2bSXiang Huang 
3528ac80d48SXiang Huang /*@C
3538ac80d48SXiang Huang   TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm
3548ac80d48SXiang Huang 
3558ac80d48SXiang Huang   Collective on Tao
3568ac80d48SXiang Huang 
3578ac80d48SXiang Huang   Level: developer
3588ac80d48SXiang Huang 
3598ac80d48SXiang Huang   Input Parameters:
3608ac80d48SXiang Huang +  tao - the Tao solver context
3618ac80d48SXiang Huang -  epsilon - L1-norm smooth approximation parameter
3628ac80d48SXiang Huang @*/
3638ac80d48SXiang Huang PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao, PetscReal epsilon)
3648ac80d48SXiang Huang {
3658ac80d48SXiang Huang   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
3668ac80d48SXiang Huang 
3678ac80d48SXiang Huang   /* Initialize epsilon here */
3688ac80d48SXiang Huang 
3698ac80d48SXiang Huang   PetscFunctionBegin;
3708ac80d48SXiang Huang   gn->epsilon = epsilon;
3718ac80d48SXiang Huang   PetscFunctionReturn(0);
3728ac80d48SXiang Huang }
373*8e85b1b3SXiang Huang 
374*8e85b1b3SXiang Huang /*@C
375*8e85b1b3SXiang Huang    TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem)
376*8e85b1b3SXiang Huang 
377*8e85b1b3SXiang Huang    Input Parameters:
378*8e85b1b3SXiang Huang +  tao  - the Tao context
379*8e85b1b3SXiang Huang .  dict - the user specified dictionary matrix
380*8e85b1b3SXiang Huang 
381*8e85b1b3SXiang Huang     Level: developer
382*8e85b1b3SXiang Huang @*/
383*8e85b1b3SXiang Huang PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao, Mat dict)
384*8e85b1b3SXiang Huang {
385*8e85b1b3SXiang Huang   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
386*8e85b1b3SXiang Huang   PetscErrorCode ierr;
387*8e85b1b3SXiang Huang   PetscFunctionBegin;
388*8e85b1b3SXiang Huang   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
389*8e85b1b3SXiang Huang   if (dict) {
390*8e85b1b3SXiang Huang     PetscValidHeaderSpecific(dict,MAT_CLASSID,2);
391*8e85b1b3SXiang Huang     /*PetscCheckSameComm(tao,1,dict,2);*/
392*8e85b1b3SXiang Huang     ierr = PetscObjectReference((PetscObject)dict);CHKERRQ(ierr);
393*8e85b1b3SXiang Huang   }
394*8e85b1b3SXiang Huang   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
395*8e85b1b3SXiang Huang   gn->D = dict;  /* We allow to set a null dictionary, which means we just use default identity matrix? */
396*8e85b1b3SXiang Huang   PetscFunctionReturn(0);
397*8e85b1b3SXiang Huang }
398*8e85b1b3SXiang Huang 
399*8e85b1b3SXiang Huang /* XH:
400*8e85b1b3SXiang Huang Changed TaoBRGNSetTikhonovLambda --> TaoBRGNSetL1RegularizerWeight  in brgn.c, peststao.h, and zbrgnf.c.
401*8e85b1b3SXiang Huang Added TaoBRGNSetL1SmoothEpsilon by following TaoBRGNSetL1RegularizerWeight.
402*8e85b1b3SXiang Huang Added TaoBRGNSetDictionaryMatrix by following TaoBRGNSetL1RegularizerWeight
4037cea06e1SXiang Huang  Maybe change D*x to D(x), and  A*x to A(x) as function handle
4047cea06e1SXiang Huang  Maybe need to also keep y = D*x, to avoid duplicate frequent computation of D*x
4057cea06e1SXiang Huang  */