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