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