xref: /petsc/src/tao/leastsquares/impls/brgn/brgn.c (revision e5a36eccef3d6b83a2c625c30d0dfd5adb4001f2)
1 #include <../src/tao/leastsquares/impls/brgn/brgn.h> /*I "petsctao.h" I*/
2 
3 #define BRGN_REGULARIZATION_USER    0
4 #define BRGN_REGULARIZATION_L2PROX  1
5 #define BRGN_REGULARIZATION_L1DICT  2
6 #define BRGN_REGULARIZATION_TYPES   3
7 
8 static const char *BRGN_REGULARIZATION_TABLE[64] = {"user","l2prox","l1dict"};
9 
10 static PetscErrorCode GNHessianProd(Mat H,Vec in,Vec out)
11 {
12   TAO_BRGN              *gn;
13   PetscErrorCode        ierr;
14 
15   PetscFunctionBegin;
16   ierr = MatShellGetContext(H,&gn);CHKERRQ(ierr);
17   ierr = MatMult(gn->subsolver->ls_jac,in,gn->r_work);CHKERRQ(ierr);
18   ierr = MatMultTranspose(gn->subsolver->ls_jac,gn->r_work,out);CHKERRQ(ierr);
19   switch (gn->reg_type) {
20   case BRGN_REGULARIZATION_USER:
21     ierr = MatMult(gn->Hreg,in,gn->x_work);CHKERRQ(ierr);
22     ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr);
23     break;
24   case BRGN_REGULARIZATION_L2PROX:
25     ierr = VecAXPY(out,gn->lambda,in);CHKERRQ(ierr);
26     break;
27   case BRGN_REGULARIZATION_L1DICT:
28     /* out = out + lambda*D'*(diag.*(D*in)) */
29     if (gn->D) {
30       ierr = MatMult(gn->D,in,gn->y);CHKERRQ(ierr);/* y = D*in */
31     } else {
32       ierr = VecCopy(in,gn->y);CHKERRQ(ierr);
33     }
34     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 */
35     if (gn->D) {
36       ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr); /* x_work = D'*(diag.*(D*in)) */
37     } else {
38       ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr);
39     }
40     ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr);
41     break;
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 static PetscErrorCode GNObjectiveGradientEval(Tao tao,Vec X,PetscReal *fcn,Vec G,void *ptr)
47 {
48   TAO_BRGN              *gn = (TAO_BRGN *)ptr;
49   PetscInt              K;                    /* dimension of D*X */
50   PetscScalar           yESum;
51   PetscErrorCode        ierr;
52   PetscReal             f_reg;
53 
54   PetscFunctionBegin;
55   /* compute objective *fcn*/
56   /* compute first term 0.5*||ls_res||_2^2 */
57   ierr = TaoComputeResidual(tao,X,tao->ls_res);CHKERRQ(ierr);
58   ierr = VecDot(tao->ls_res,tao->ls_res,fcn);CHKERRQ(ierr);
59   *fcn *= 0.5;
60   /* compute gradient G */
61   ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr);
62   ierr = MatMultTranspose(tao->ls_jac,tao->ls_res,G);CHKERRQ(ierr);
63   /* add the regularization contribution */
64   switch (gn->reg_type) {
65   case BRGN_REGULARIZATION_USER:
66     ierr = (*gn->regularizerobjandgrad)(tao,X,&f_reg,gn->x_work,gn->reg_obj_ctx);CHKERRQ(ierr);
67     *fcn += gn->lambda*f_reg;
68     ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr);
69     break;
70   case BRGN_REGULARIZATION_L2PROX:
71     /* compute f = f + lambda*0.5*(xk - xkm1)'*(xk - xkm1) */
72     ierr = VecAXPBYPCZ(gn->x_work,1.0,-1.0,0.0,X,gn->x_old);CHKERRQ(ierr);
73     ierr = VecDot(gn->x_work,gn->x_work,&f_reg);CHKERRQ(ierr);
74     *fcn += gn->lambda*0.5*f_reg;
75     /* compute G = G + lambda*(xk - xkm1) */
76     ierr = VecAXPBYPCZ(G,gn->lambda,-gn->lambda,1.0,X,gn->x_old);CHKERRQ(ierr);
77     break;
78   case BRGN_REGULARIZATION_L1DICT:
79     /* compute f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x*/
80     if (gn->D) {
81       ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */
82     } else {
83       ierr = VecCopy(X,gn->y);CHKERRQ(ierr);
84     }
85     ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr);
86     ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
87     ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr);  /* gn->y_work = sqrt(y.^2+epsilon^2) */
88     ierr = VecSum(gn->y_work,&yESum);CHKERRQ(ierr);CHKERRQ(ierr);
89     ierr = VecGetSize(gn->y,&K);CHKERRQ(ierr);
90     *fcn += gn->lambda*(yESum - K*gn->epsilon);
91     /* compute G = G + lambda*D'*(y./sqrt(y.^2+epsilon^2)),where y = D*x */
92     ierr = VecPointwiseDivide(gn->y_work,gn->y,gn->y_work);CHKERRQ(ierr); /* reuse y_work = y./sqrt(y.^2+epsilon^2) */
93     if (gn->D) {
94       ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr);
95     } else {
96       ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr);
97     }
98     ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr);
99     break;
100   }
101   PetscFunctionReturn(0);
102 }
103 
104 static PetscErrorCode GNComputeHessian(Tao tao,Vec X,Mat H,Mat Hpre,void *ptr)
105 {
106   TAO_BRGN              *gn = (TAO_BRGN *)ptr;
107   PetscErrorCode ierr;
108 
109   PetscFunctionBegin;
110   ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr);
111 
112   switch (gn->reg_type) {
113   case BRGN_REGULARIZATION_USER:
114     ierr = (*gn->regularizerhessian)(tao,X,gn->Hreg,gn->reg_hess_ctx);CHKERRQ(ierr);
115     break;
116   case BRGN_REGULARIZATION_L2PROX:
117     break;
118   case BRGN_REGULARIZATION_L1DICT:
119     /* 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 */
120     if (gn->D) {
121       ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */
122     } else {
123       ierr = VecCopy(X,gn->y);CHKERRQ(ierr);
124     }
125     ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr);
126     ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
127     ierr = VecCopy(gn->y_work,gn->diag);CHKERRQ(ierr);                  /* gn->diag = y.^2+epsilon^2 */
128     ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr);                        /* gn->y_work = sqrt(y.^2+epsilon^2) */
129     ierr = VecPointwiseMult(gn->diag,gn->y_work,gn->diag);CHKERRQ(ierr);/* gn->diag = sqrt(y.^2+epsilon^2).^3 */
130     ierr = VecReciprocal(gn->diag);CHKERRQ(ierr);
131     ierr = VecScale(gn->diag,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
132     break;
133   }
134   PetscFunctionReturn(0);
135 }
136 
137 static PetscErrorCode GNHookFunction(Tao tao,PetscInt iter, void *ctx)
138 {
139   TAO_BRGN              *gn = (TAO_BRGN *)ctx;
140   PetscErrorCode        ierr;
141 
142   PetscFunctionBegin;
143   /* Update basic tao information from the subsolver */
144   gn->parent->nfuncs = tao->nfuncs;
145   gn->parent->ngrads = tao->ngrads;
146   gn->parent->nfuncgrads = tao->nfuncgrads;
147   gn->parent->nhess = tao->nhess;
148   gn->parent->niter = tao->niter;
149   gn->parent->ksp_its = tao->ksp_its;
150   gn->parent->ksp_tot_its = tao->ksp_tot_its;
151   ierr = TaoGetConvergedReason(tao,&gn->parent->reason);CHKERRQ(ierr);
152   /* Update the solution vectors */
153   if (iter == 0) {
154     ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr);
155   } else {
156     ierr = VecCopy(tao->solution,gn->x_old);CHKERRQ(ierr);
157     ierr = VecCopy(tao->solution,gn->parent->solution);CHKERRQ(ierr);
158   }
159   /* Update the gradient */
160   ierr = VecCopy(tao->gradient,gn->parent->gradient);CHKERRQ(ierr);
161   /* Call general purpose update function */
162   if (gn->parent->ops->update) {
163     ierr = (*gn->parent->ops->update)(gn->parent,gn->parent->niter,gn->parent->user_update);CHKERRQ(ierr);
164   }
165   PetscFunctionReturn(0);
166 }
167 
168 static PetscErrorCode TaoSolve_BRGN(Tao tao)
169 {
170   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
171   PetscErrorCode        ierr;
172 
173   PetscFunctionBegin;
174   ierr = TaoSolve(gn->subsolver);CHKERRQ(ierr);
175   /* Update basic tao information from the subsolver */
176   tao->nfuncs = gn->subsolver->nfuncs;
177   tao->ngrads = gn->subsolver->ngrads;
178   tao->nfuncgrads = gn->subsolver->nfuncgrads;
179   tao->nhess = gn->subsolver->nhess;
180   tao->niter = gn->subsolver->niter;
181   tao->ksp_its = gn->subsolver->ksp_its;
182   tao->ksp_tot_its = gn->subsolver->ksp_tot_its;
183   ierr = TaoGetConvergedReason(gn->subsolver,&tao->reason);CHKERRQ(ierr);
184   /* Update vectors */
185   ierr = VecCopy(gn->subsolver->solution,tao->solution);CHKERRQ(ierr);
186   ierr = VecCopy(gn->subsolver->gradient,tao->gradient);CHKERRQ(ierr);
187   PetscFunctionReturn(0);
188 }
189 
190 static PetscErrorCode TaoSetFromOptions_BRGN(PetscOptionItems *PetscOptionsObject,Tao tao)
191 {
192   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
193   PetscErrorCode        ierr;
194 
195   PetscFunctionBegin;
196   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);
197   ierr = PetscOptionsReal("-tao_brgn_regularizer_weight","regularizer weight (default 1e-4)","",gn->lambda,&gn->lambda,NULL);CHKERRQ(ierr);
198   ierr = PetscOptionsReal("-tao_brgn_l1_smooth_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);
199   ierr = PetscOptionsEList("-tao_brgn_regularization_type","regularization type", "",BRGN_REGULARIZATION_TABLE,BRGN_REGULARIZATION_TYPES,BRGN_REGULARIZATION_TABLE[gn->reg_type],&gn->reg_type,NULL);CHKERRQ(ierr);
200   ierr = PetscOptionsTail();CHKERRQ(ierr);
201   ierr = TaoSetFromOptions(gn->subsolver);CHKERRQ(ierr);
202   PetscFunctionReturn(0);
203 }
204 
205 static PetscErrorCode TaoView_BRGN(Tao tao,PetscViewer viewer)
206 {
207   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
208   PetscErrorCode        ierr;
209 
210   PetscFunctionBegin;
211   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
212   ierr = TaoView(gn->subsolver,viewer);CHKERRQ(ierr);
213   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
214   PetscFunctionReturn(0);
215 }
216 
217 static PetscErrorCode TaoSetUp_BRGN(Tao tao)
218 {
219   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
220   PetscErrorCode        ierr;
221   PetscBool             is_bnls,is_bntr,is_bntl;
222   PetscInt              i,n,N,K; /* dict has size K*N*/
223 
224   PetscFunctionBegin;
225   if (!tao->ls_res) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualRoutine() must be called before setup!");
226   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNLS,&is_bnls);CHKERRQ(ierr);
227   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTR,&is_bntr);CHKERRQ(ierr);
228   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTL,&is_bntl);CHKERRQ(ierr);
229   if ((is_bnls || is_bntr || is_bntl) && !tao->ls_jac) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualJacobianRoutine() must be called before setup!");
230   if (!tao->gradient) {
231     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
232   }
233   if (!gn->x_work) {
234     ierr = VecDuplicate(tao->solution,&gn->x_work);CHKERRQ(ierr);
235   }
236   if (!gn->r_work) {
237     ierr = VecDuplicate(tao->ls_res,&gn->r_work);CHKERRQ(ierr);
238   }
239   if (!gn->x_old) {
240     ierr = VecDuplicate(tao->solution,&gn->x_old);CHKERRQ(ierr);
241     ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr);
242   }
243 
244   if (BRGN_REGULARIZATION_L1DICT == gn->reg_type) {
245     if (gn->D) {
246       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 */
247     } else {
248       ierr = VecGetSize(tao->solution,&K);CHKERRQ(ierr); /* If user does not setup dict matrix, use identiy matrix, K=N */
249     }
250     if (!gn->y) {
251       ierr = VecCreate(PETSC_COMM_SELF,&gn->y);CHKERRQ(ierr);
252       ierr = VecSetSizes(gn->y,PETSC_DECIDE,K);CHKERRQ(ierr);
253       ierr = VecSetFromOptions(gn->y);CHKERRQ(ierr);
254       ierr = VecSet(gn->y,0.0);CHKERRQ(ierr);
255 
256     }
257     if (!gn->y_work) {
258       ierr = VecDuplicate(gn->y,&gn->y_work);CHKERRQ(ierr);
259     }
260     if (!gn->diag) {
261       ierr = VecDuplicate(gn->y,&gn->diag);CHKERRQ(ierr);
262       ierr = VecSet(gn->diag,0.0);CHKERRQ(ierr);
263     }
264   }
265 
266   if (!tao->setupcalled) {
267     /* Hessian setup */
268     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
269     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
270     ierr = MatSetSizes(gn->H,n,n,N,N);CHKERRQ(ierr);
271     ierr = MatSetType(gn->H,MATSHELL);CHKERRQ(ierr);
272     ierr = MatSetUp(gn->H);CHKERRQ(ierr);
273     ierr = MatShellSetOperation(gn->H,MATOP_MULT,(void (*)(void))GNHessianProd);CHKERRQ(ierr);
274     ierr = MatShellSetContext(gn->H,(void*)gn);CHKERRQ(ierr);
275     /* Subsolver setup,include initial vector and dicttionary D */
276     ierr = TaoSetUpdate(gn->subsolver,GNHookFunction,(void*)gn);CHKERRQ(ierr);
277     ierr = TaoSetInitialVector(gn->subsolver,tao->solution);CHKERRQ(ierr);
278     if (tao->bounded) {
279       ierr = TaoSetVariableBounds(gn->subsolver,tao->XL,tao->XU);CHKERRQ(ierr);
280     }
281     ierr = TaoSetResidualRoutine(gn->subsolver,tao->ls_res,tao->ops->computeresidual,tao->user_lsresP);CHKERRQ(ierr);
282     ierr = TaoSetJacobianResidualRoutine(gn->subsolver,tao->ls_jac,tao->ls_jac,tao->ops->computeresidualjacobian,tao->user_lsjacP);CHKERRQ(ierr);
283     ierr = TaoSetObjectiveAndGradientRoutine(gn->subsolver,GNObjectiveGradientEval,(void*)gn);CHKERRQ(ierr);
284     ierr = TaoSetHessianRoutine(gn->subsolver,gn->H,gn->H,GNComputeHessian,(void*)gn);CHKERRQ(ierr);
285     /* Propagate some options down */
286     ierr = TaoSetTolerances(gn->subsolver,tao->gatol,tao->grtol,tao->gttol);CHKERRQ(ierr);
287     ierr = TaoSetMaximumIterations(gn->subsolver,tao->max_it);CHKERRQ(ierr);
288     ierr = TaoSetMaximumFunctionEvaluations(gn->subsolver,tao->max_funcs);CHKERRQ(ierr);
289     for (i=0; i<tao->numbermonitors; ++i) {
290       ierr = TaoSetMonitor(gn->subsolver,tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i]);CHKERRQ(ierr);
291       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
292     }
293     ierr = TaoSetUp(gn->subsolver);CHKERRQ(ierr);
294   }
295   PetscFunctionReturn(0);
296 }
297 
298 static PetscErrorCode TaoDestroy_BRGN(Tao tao)
299 {
300   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
301   PetscErrorCode        ierr;
302 
303   PetscFunctionBegin;
304   if (tao->setupcalled) {
305     ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr);
306     ierr = VecDestroy(&gn->x_work);CHKERRQ(ierr);
307     ierr = VecDestroy(&gn->r_work);CHKERRQ(ierr);
308     ierr = VecDestroy(&gn->x_old);CHKERRQ(ierr);
309     ierr = VecDestroy(&gn->diag);CHKERRQ(ierr);
310     ierr = VecDestroy(&gn->y);CHKERRQ(ierr);
311     ierr = VecDestroy(&gn->y_work);CHKERRQ(ierr);
312   }
313   ierr = MatDestroy(&gn->H);CHKERRQ(ierr);
314   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
315   ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr);
316   ierr = TaoDestroy(&gn->subsolver);CHKERRQ(ierr);
317   gn->parent = NULL;
318   ierr = PetscFree(tao->data);CHKERRQ(ierr);
319   PetscFunctionReturn(0);
320 }
321 
322 /*MC
323   TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares
324             problems with bound constraints. This algorithm is a thin wrapper around TAOBNTL
325             that constructs the Gauss-Newton problem with the user-provided least-squares
326             residual and Jacobian. The algorithm offers both an L2-norm proximal point ("l2prox")
327             regularizer, and a L1-norm dictionary regularizer ("l1dict"), where we approximate the
328             L1-norm ||x||_1 by sum_i(sqrt(x_i^2+epsilon^2)-epsilon) with a small positive number epsilon.
329             The user can also provide own regularization function.
330 
331   Options Database Keys:
332   + -tao_brgn_regularization_type - regularization type ("user", "l2prox", "l1dict") (default "l2prox")
333   . -tao_brgn_regularizer_weight  - regularizer weight (default 1e-4)
334   - -tao_brgn_l1_smooth_epsilon   - L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6)
335 
336   Level: beginner
337 M*/
338 PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao)
339 {
340   TAO_BRGN       *gn;
341   PetscErrorCode ierr;
342 
343   PetscFunctionBegin;
344   ierr = PetscNewLog(tao,&gn);CHKERRQ(ierr);
345 
346   tao->ops->destroy = TaoDestroy_BRGN;
347   tao->ops->setup = TaoSetUp_BRGN;
348   tao->ops->setfromoptions = TaoSetFromOptions_BRGN;
349   tao->ops->view = TaoView_BRGN;
350   tao->ops->solve = TaoSolve_BRGN;
351 
352   tao->data = (void*)gn;
353   gn->reg_type = BRGN_REGULARIZATION_L2PROX;
354   gn->lambda = 1e-4;
355   gn->epsilon = 1e-6;
356   gn->parent = tao;
357 
358   ierr = MatCreate(PetscObjectComm((PetscObject)tao),&gn->H);CHKERRQ(ierr);
359   ierr = MatSetOptionsPrefix(gn->H,"tao_brgn_hessian_");CHKERRQ(ierr);
360 
361   ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&gn->subsolver);CHKERRQ(ierr);
362   ierr = TaoSetType(gn->subsolver,TAOBNLS);CHKERRQ(ierr);
363   ierr = TaoSetOptionsPrefix(gn->subsolver,"tao_brgn_subsolver_");CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 }
366 
367 /*@
368   TaoBRGNGetSubsolver - Get the pointer to the subsolver inside BRGN
369 
370   Collective on Tao
371 
372   Level: advanced
373 
374   Input Parameters:
375 +  tao - the Tao solver context
376 -  subsolver - the Tao sub-solver context
377 @*/
378 PetscErrorCode TaoBRGNGetSubsolver(Tao tao,Tao *subsolver)
379 {
380   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
381 
382   PetscFunctionBegin;
383   *subsolver = gn->subsolver;
384   PetscFunctionReturn(0);
385 }
386 
387 /*@
388   TaoBRGNSetRegularizerWeight - Set the regularizer weight for the Gauss-Newton least-squares algorithm
389 
390   Collective on Tao
391 
392   Input Parameters:
393 +  tao - the Tao solver context
394 -  lambda - L1-norm regularizer weight
395 
396   Level: beginner
397 @*/
398 PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao,PetscReal lambda)
399 {
400   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
401 
402   /* Initialize lambda here */
403 
404   PetscFunctionBegin;
405   gn->lambda = lambda;
406   PetscFunctionReturn(0);
407 }
408 
409 /*@
410   TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm
411 
412   Collective on Tao
413 
414   Input Parameters:
415 +  tao - the Tao solver context
416 -  epsilon - L1-norm smooth approximation parameter
417 
418   Level: advanced
419 @*/
420 PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao,PetscReal epsilon)
421 {
422   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
423 
424   /* Initialize epsilon here */
425 
426   PetscFunctionBegin;
427   gn->epsilon = epsilon;
428   PetscFunctionReturn(0);
429 }
430 
431 /*@
432    TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem)
433 
434    Input Parameters:
435 +  tao  - the Tao context
436 .  dict - the user specified dictionary matrix.  We allow to set a null dictionary, which means identity matrix by default
437 
438     Level: advanced
439 @*/
440 PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao,Mat dict)
441 {
442   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
443   PetscErrorCode ierr;
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
446   if (dict) {
447     PetscValidHeaderSpecific(dict,MAT_CLASSID,2);
448     PetscCheckSameComm(tao,1,dict,2);
449     ierr = PetscObjectReference((PetscObject)dict);CHKERRQ(ierr);
450   }
451   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
452   gn->D = dict;
453   PetscFunctionReturn(0);
454 }
455 
456 /*@C
457    TaoBRGNSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back
458    function into the algorithm.
459 
460    Input Parameters:
461    + tao - the Tao context
462    . func - function pointer for the regularizer value and gradient evaluation
463    - ctx - user context for the regularizer
464 
465    Level: advanced
466 @*/
467 PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (*func)(Tao,Vec,PetscReal *,Vec,void*),void *ctx)
468 {
469   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
470 
471   PetscFunctionBegin;
472   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
473   if (ctx) {
474     gn->reg_obj_ctx = ctx;
475   }
476   if (func) {
477     gn->regularizerobjandgrad = func;
478   }
479   PetscFunctionReturn(0);
480 }
481 
482 /*@C
483    TaoBRGNSetRegularizerHessianRoutine - Sets the user-defined regularizer call-back
484    function into the algorithm.
485 
486    Input Parameters:
487    + tao - the Tao context
488    . Hreg - user-created matrix for the Hessian of the regularization term
489    . func - function pointer for the regularizer Hessian evaluation
490    - ctx - user context for the regularizer Hessian
491 
492    Level: advanced
493 @*/
494 PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao,Mat Hreg,PetscErrorCode (*func)(Tao,Vec,Mat,void*),void *ctx)
495 {
496   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
497   PetscErrorCode ierr;
498 
499   PetscFunctionBegin;
500   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
501   if (Hreg) {
502     PetscValidHeaderSpecific(Hreg,MAT_CLASSID,2);
503     PetscCheckSameComm(tao,1,Hreg,2);
504   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONG,"NULL Hessian detected! User must provide valid Hessian for the regularizer.");
505   if (ctx) {
506     gn->reg_hess_ctx = ctx;
507   }
508   if (func) {
509     gn->regularizerhessian = func;
510   }
511   if (Hreg) {
512     ierr = PetscObjectReference((PetscObject)Hreg);CHKERRQ(ierr);
513     ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr);
514     gn->Hreg = Hreg;
515   }
516   PetscFunctionReturn(0);
517 }
518