xref: /petsc/src/tao/leastsquares/impls/brgn/brgn.c (revision f2ed2dc71a2ab9ffda85eae8afa0cbea9ed570de)
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 
160   switch (gn->reg_type) {
161   case BRGN_REGULARIZATION_USER:
162     ierr = (*gn->regularizerhessian)(tao,X,gn->Hreg,gn->reg_hess_ctx);CHKERRQ(ierr);
163     break;
164   case BRGN_REGULARIZATION_L2PURE:
165     break;
166   case BRGN_REGULARIZATION_L2PROX:
167     break;
168   case BRGN_REGULARIZATION_L1DICT:
169     /* 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 */
170     if (gn->D) {
171       ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */
172     } else {
173       ierr = VecCopy(X,gn->y);CHKERRQ(ierr);
174     }
175     ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr);
176     ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
177     ierr = VecCopy(gn->y_work,gn->diag);CHKERRQ(ierr);                  /* gn->diag = y.^2+epsilon^2 */
178     ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr);                        /* gn->y_work = sqrt(y.^2+epsilon^2) */
179     ierr = VecPointwiseMult(gn->diag,gn->y_work,gn->diag);CHKERRQ(ierr);/* gn->diag = sqrt(y.^2+epsilon^2).^3 */
180     ierr = VecReciprocal(gn->diag);CHKERRQ(ierr);
181     ierr = VecScale(gn->diag,gn->epsilon*gn->epsilon);CHKERRQ(ierr);
182     break;
183   case BRGN_REGULARIZATION_LM:
184     /* compute diagonal of J^T J */
185     ierr = MatGetSize(gn->parent->ls_jac,NULL,&n);CHKERRQ(ierr);
186     ierr = PetscMalloc1(n,&cnorms);CHKERRQ(ierr);
187     ierr = MatGetColumnNorms(gn->parent->ls_jac,NORM_2,cnorms);CHKERRQ(ierr);
188     ierr = MatGetOwnershipRangeColumn(gn->parent->ls_jac,&cstart,&cend);CHKERRQ(ierr);
189     ierr = VecGetArray(gn->diag,&diag_ary);CHKERRQ(ierr);
190     for (i = 0; i < cend-cstart; i++) {
191       diag_ary[i] = cnorms[cstart+i] * cnorms[cstart+i];
192     }
193     ierr = VecRestoreArray(gn->diag,&diag_ary);CHKERRQ(ierr);
194     ierr = PetscFree(cnorms);CHKERRQ(ierr);
195     ierr = ComputeDamping(gn);CHKERRQ(ierr);
196     break;
197   }
198   PetscFunctionReturn(0);
199 }
200 
201 static PetscErrorCode GNHookFunction(Tao tao,PetscInt iter, void *ctx)
202 {
203   TAO_BRGN              *gn = (TAO_BRGN *)ctx;
204   PetscErrorCode        ierr;
205 
206   PetscFunctionBegin;
207   /* Update basic tao information from the subsolver */
208   gn->parent->nfuncs = tao->nfuncs;
209   gn->parent->ngrads = tao->ngrads;
210   gn->parent->nfuncgrads = tao->nfuncgrads;
211   gn->parent->nhess = tao->nhess;
212   gn->parent->niter = tao->niter;
213   gn->parent->ksp_its = tao->ksp_its;
214   gn->parent->ksp_tot_its = tao->ksp_tot_its;
215   gn->parent->fc = tao->fc;
216   ierr = TaoGetConvergedReason(tao,&gn->parent->reason);CHKERRQ(ierr);
217   /* Update the solution vectors */
218   if (iter == 0) {
219     ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr);
220   } else {
221     ierr = VecCopy(tao->solution,gn->x_old);CHKERRQ(ierr);
222     ierr = VecCopy(tao->solution,gn->parent->solution);CHKERRQ(ierr);
223   }
224   /* Update the gradient */
225   ierr = VecCopy(tao->gradient,gn->parent->gradient);CHKERRQ(ierr);
226 
227   /* Update damping parameter for LM */
228   if (gn->reg_type == BRGN_REGULARIZATION_LM) {
229     if (iter > 0) {
230       if (gn->fc_old > tao->fc) {
231         gn->lambda = gn->lambda * gn->downhill_lambda_change;
232       } else {
233         /* uphill step */
234         gn->lambda = gn->lambda * gn->uphill_lambda_change;
235       }
236     }
237     gn->fc_old = tao->fc;
238   }
239 
240   /* Call general purpose update function */
241   if (gn->parent->ops->update) {
242     ierr = (*gn->parent->ops->update)(gn->parent,gn->parent->niter,gn->parent->user_update);CHKERRQ(ierr);
243   }
244   PetscFunctionReturn(0);
245 }
246 
247 static PetscErrorCode TaoSolve_BRGN(Tao tao)
248 {
249   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
250   PetscErrorCode        ierr;
251 
252   PetscFunctionBegin;
253   ierr = TaoSolve(gn->subsolver);CHKERRQ(ierr);
254   /* Update basic tao information from the subsolver */
255   tao->nfuncs = gn->subsolver->nfuncs;
256   tao->ngrads = gn->subsolver->ngrads;
257   tao->nfuncgrads = gn->subsolver->nfuncgrads;
258   tao->nhess = gn->subsolver->nhess;
259   tao->niter = gn->subsolver->niter;
260   tao->ksp_its = gn->subsolver->ksp_its;
261   tao->ksp_tot_its = gn->subsolver->ksp_tot_its;
262   ierr = TaoGetConvergedReason(gn->subsolver,&tao->reason);CHKERRQ(ierr);
263   /* Update vectors */
264   ierr = VecCopy(gn->subsolver->solution,tao->solution);CHKERRQ(ierr);
265   ierr = VecCopy(gn->subsolver->gradient,tao->gradient);CHKERRQ(ierr);
266   PetscFunctionReturn(0);
267 }
268 
269 static PetscErrorCode TaoSetFromOptions_BRGN(PetscOptionItems *PetscOptionsObject,Tao tao)
270 {
271   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
272   TaoLineSearch         ls;
273   PetscErrorCode        ierr;
274 
275   PetscFunctionBegin;
276   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);
277   ierr = PetscOptionsReal("-tao_brgn_regularizer_weight","regularizer weight (default 1e-4)","",gn->lambda,&gn->lambda,NULL);CHKERRQ(ierr);
278   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);
279   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);
280   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);
281   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);
282   ierr = PetscOptionsTail();CHKERRQ(ierr);
283   /* set unit line search direction as the default when using the lm regularizer */
284   if (gn->reg_type == BRGN_REGULARIZATION_LM) {
285     ierr = TaoGetLineSearch(gn->subsolver,&ls);CHKERRQ(ierr);
286     ierr = TaoLineSearchSetType(ls,TAOLINESEARCHUNIT);CHKERRQ(ierr);
287   }
288   ierr = TaoSetFromOptions(gn->subsolver);CHKERRQ(ierr);
289   PetscFunctionReturn(0);
290 }
291 
292 static PetscErrorCode TaoView_BRGN(Tao tao,PetscViewer viewer)
293 {
294   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
295   PetscErrorCode        ierr;
296 
297   PetscFunctionBegin;
298   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
299   ierr = TaoView(gn->subsolver,viewer);CHKERRQ(ierr);
300   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303 
304 static PetscErrorCode TaoSetUp_BRGN(Tao tao)
305 {
306   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
307   PetscErrorCode        ierr;
308   PetscBool             is_bnls,is_bntr,is_bntl;
309   PetscInt              i,n,N,K; /* dict has size K*N*/
310 
311   PetscFunctionBegin;
312   if (!tao->ls_res) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualRoutine() must be called before setup!");
313   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNLS,&is_bnls);CHKERRQ(ierr);
314   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTR,&is_bntr);CHKERRQ(ierr);
315   ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTL,&is_bntl);CHKERRQ(ierr);
316   if ((is_bnls || is_bntr || is_bntl) && !tao->ls_jac) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualJacobianRoutine() must be called before setup!");
317   if (!tao->gradient) {
318     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
319   }
320   if (!gn->x_work) {
321     ierr = VecDuplicate(tao->solution,&gn->x_work);CHKERRQ(ierr);
322   }
323   if (!gn->r_work) {
324     ierr = VecDuplicate(tao->ls_res,&gn->r_work);CHKERRQ(ierr);
325   }
326   if (!gn->x_old) {
327     ierr = VecDuplicate(tao->solution,&gn->x_old);CHKERRQ(ierr);
328     ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr);
329   }
330 
331   if (BRGN_REGULARIZATION_L1DICT == gn->reg_type) {
332     if (gn->D) {
333       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 */
334     } else {
335       ierr = VecGetSize(tao->solution,&K);CHKERRQ(ierr); /* If user does not setup dict matrix, use identiy matrix, K=N */
336     }
337     if (!gn->y) {
338       ierr = VecCreate(PETSC_COMM_SELF,&gn->y);CHKERRQ(ierr);
339       ierr = VecSetSizes(gn->y,PETSC_DECIDE,K);CHKERRQ(ierr);
340       ierr = VecSetFromOptions(gn->y);CHKERRQ(ierr);
341       ierr = VecSet(gn->y,0.0);CHKERRQ(ierr);
342 
343     }
344     if (!gn->y_work) {
345       ierr = VecDuplicate(gn->y,&gn->y_work);CHKERRQ(ierr);
346     }
347     if (!gn->diag) {
348       ierr = VecDuplicate(gn->y,&gn->diag);CHKERRQ(ierr);
349       ierr = VecSet(gn->diag,0.0);CHKERRQ(ierr);
350     }
351   }
352   if (BRGN_REGULARIZATION_LM == gn->reg_type) {
353     if (!gn->diag) {
354       ierr = MatCreateVecs(gn->parent->ls_jac,&gn->diag,NULL);CHKERRQ(ierr);
355     }
356     if (!gn->damping) {
357       ierr = MatCreateVecs(gn->parent->ls_jac,&gn->damping,NULL);CHKERRQ(ierr);
358     }
359   }
360 
361   if (!tao->setupcalled) {
362     /* Hessian setup */
363     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
364     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
365     ierr = MatSetSizes(gn->H,n,n,N,N);CHKERRQ(ierr);
366     ierr = MatSetType(gn->H,MATSHELL);CHKERRQ(ierr);
367     ierr = MatSetUp(gn->H);CHKERRQ(ierr);
368     ierr = MatShellSetOperation(gn->H,MATOP_MULT,(void (*)(void))GNHessianProd);CHKERRQ(ierr);
369     ierr = MatShellSetContext(gn->H,(void*)gn);CHKERRQ(ierr);
370     /* Subsolver setup,include initial vector and dicttionary D */
371     ierr = TaoSetUpdate(gn->subsolver,GNHookFunction,(void*)gn);CHKERRQ(ierr);
372     ierr = TaoSetInitialVector(gn->subsolver,tao->solution);CHKERRQ(ierr);
373     if (tao->bounded) {
374       ierr = TaoSetVariableBounds(gn->subsolver,tao->XL,tao->XU);CHKERRQ(ierr);
375     }
376     ierr = TaoSetResidualRoutine(gn->subsolver,tao->ls_res,tao->ops->computeresidual,tao->user_lsresP);CHKERRQ(ierr);
377     ierr = TaoSetJacobianResidualRoutine(gn->subsolver,tao->ls_jac,tao->ls_jac,tao->ops->computeresidualjacobian,tao->user_lsjacP);CHKERRQ(ierr);
378     ierr = TaoSetObjectiveAndGradientRoutine(gn->subsolver,GNObjectiveGradientEval,(void*)gn);CHKERRQ(ierr);
379     ierr = TaoSetHessianRoutine(gn->subsolver,gn->H,gn->H,GNComputeHessian,(void*)gn);CHKERRQ(ierr);
380     /* Propagate some options down */
381     ierr = TaoSetTolerances(gn->subsolver,tao->gatol,tao->grtol,tao->gttol);CHKERRQ(ierr);
382     ierr = TaoSetMaximumIterations(gn->subsolver,tao->max_it);CHKERRQ(ierr);
383     ierr = TaoSetMaximumFunctionEvaluations(gn->subsolver,tao->max_funcs);CHKERRQ(ierr);
384     for (i=0; i<tao->numbermonitors; ++i) {
385       ierr = TaoSetMonitor(gn->subsolver,tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i]);CHKERRQ(ierr);
386       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
387     }
388     ierr = TaoSetUp(gn->subsolver);CHKERRQ(ierr);
389   }
390   PetscFunctionReturn(0);
391 }
392 
393 static PetscErrorCode TaoDestroy_BRGN(Tao tao)
394 {
395   TAO_BRGN              *gn = (TAO_BRGN *)tao->data;
396   PetscErrorCode        ierr;
397 
398   PetscFunctionBegin;
399   if (tao->setupcalled) {
400     ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr);
401     ierr = VecDestroy(&gn->x_work);CHKERRQ(ierr);
402     ierr = VecDestroy(&gn->r_work);CHKERRQ(ierr);
403     ierr = VecDestroy(&gn->x_old);CHKERRQ(ierr);
404     ierr = VecDestroy(&gn->diag);CHKERRQ(ierr);
405     ierr = VecDestroy(&gn->y);CHKERRQ(ierr);
406     ierr = VecDestroy(&gn->y_work);CHKERRQ(ierr);
407   }
408   ierr = VecDestroy(&gn->damping);CHKERRQ(ierr);
409   ierr = VecDestroy(&gn->diag);CHKERRQ(ierr);
410   ierr = MatDestroy(&gn->H);CHKERRQ(ierr);
411   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
412   ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr);
413   ierr = TaoDestroy(&gn->subsolver);CHKERRQ(ierr);
414   gn->parent = NULL;
415   ierr = PetscFree(tao->data);CHKERRQ(ierr);
416   PetscFunctionReturn(0);
417 }
418 
419 /*MC
420   TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares
421             problems with bound constraints. This algorithm is a thin wrapper around TAOBNTL
422             that constructs the Gauss-Newton problem with the user-provided least-squares
423             residual and Jacobian. The algorithm offers an L2-norm ("l2pure"), L2-norm proximal point ("l2prox")
424             regularizer, and L1-norm dictionary regularizer ("l1dict"), where we approximate the
425             L1-norm ||x||_1 by sum_i(sqrt(x_i^2+epsilon^2)-epsilon) with a small positive number epsilon.
426             Also offered is the "lm" regularizer which uses a scaled diagonal of J^T J.
427             With the "lm" regularizer, BRGN is a Levenberg-Marquardt optimizer.
428             The user can also provide own regularization function.
429 
430   Options Database Keys:
431 + -tao_brgn_regularization_type - regularization type ("user", "l2prox", "l2pure", "l1dict", "lm") (default "l2prox")
432 . -tao_brgn_regularizer_weight  - regularizer weight (default 1e-4)
433 - -tao_brgn_l1_smooth_epsilon   - L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6)
434 
435   Level: beginner
436 M*/
437 PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao)
438 {
439   TAO_BRGN       *gn;
440   PetscErrorCode ierr;
441 
442   PetscFunctionBegin;
443   ierr = PetscNewLog(tao,&gn);CHKERRQ(ierr);
444 
445   tao->ops->destroy = TaoDestroy_BRGN;
446   tao->ops->setup = TaoSetUp_BRGN;
447   tao->ops->setfromoptions = TaoSetFromOptions_BRGN;
448   tao->ops->view = TaoView_BRGN;
449   tao->ops->solve = TaoSolve_BRGN;
450 
451   tao->data = (void*)gn;
452   gn->reg_type = BRGN_REGULARIZATION_L2PROX;
453   gn->lambda = 1e-4;
454   gn->epsilon = 1e-6;
455   gn->downhill_lambda_change = 1./5.;
456   gn->uphill_lambda_change = 1.5;
457   gn->parent = tao;
458 
459   ierr = MatCreate(PetscObjectComm((PetscObject)tao),&gn->H);CHKERRQ(ierr);
460   ierr = MatSetOptionsPrefix(gn->H,"tao_brgn_hessian_");CHKERRQ(ierr);
461 
462   ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&gn->subsolver);CHKERRQ(ierr);
463   ierr = TaoSetType(gn->subsolver,TAOBNLS);CHKERRQ(ierr);
464   ierr = TaoSetOptionsPrefix(gn->subsolver,"tao_brgn_subsolver_");CHKERRQ(ierr);
465   PetscFunctionReturn(0);
466 }
467 
468 /*@
469   TaoBRGNGetSubsolver - Get the pointer to the subsolver inside BRGN
470 
471   Collective on Tao
472 
473   Level: advanced
474 
475   Input Parameters:
476 +  tao - the Tao solver context
477 -  subsolver - the Tao sub-solver context
478 @*/
479 PetscErrorCode TaoBRGNGetSubsolver(Tao tao,Tao *subsolver)
480 {
481   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
482 
483   PetscFunctionBegin;
484   *subsolver = gn->subsolver;
485   PetscFunctionReturn(0);
486 }
487 
488 /*@
489   TaoBRGNSetRegularizerWeight - Set the regularizer weight for the Gauss-Newton least-squares algorithm
490 
491   Collective on Tao
492 
493   Input Parameters:
494 +  tao - the Tao solver context
495 -  lambda - L1-norm regularizer weight
496 
497   Level: beginner
498 @*/
499 PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao,PetscReal lambda)
500 {
501   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
502 
503   /* Initialize lambda here */
504 
505   PetscFunctionBegin;
506   gn->lambda = lambda;
507   PetscFunctionReturn(0);
508 }
509 
510 /*@
511   TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm
512 
513   Collective on Tao
514 
515   Input Parameters:
516 +  tao - the Tao solver context
517 -  epsilon - L1-norm smooth approximation parameter
518 
519   Level: advanced
520 @*/
521 PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao,PetscReal epsilon)
522 {
523   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
524 
525   /* Initialize epsilon here */
526 
527   PetscFunctionBegin;
528   gn->epsilon = epsilon;
529   PetscFunctionReturn(0);
530 }
531 
532 /*@
533    TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem)
534 
535    Input Parameters:
536 +  tao  - the Tao context
537 -  dict - the user specified dictionary matrix.  We allow to set a null dictionary, which means identity matrix by default
538 
539     Level: advanced
540 @*/
541 PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao,Mat dict)
542 {
543   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
544   PetscErrorCode ierr;
545   PetscFunctionBegin;
546   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
547   if (dict) {
548     PetscValidHeaderSpecific(dict,MAT_CLASSID,2);
549     PetscCheckSameComm(tao,1,dict,2);
550     ierr = PetscObjectReference((PetscObject)dict);CHKERRQ(ierr);
551   }
552   ierr = MatDestroy(&gn->D);CHKERRQ(ierr);
553   gn->D = dict;
554   PetscFunctionReturn(0);
555 }
556 
557 /*@C
558    TaoBRGNSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back
559    function into the algorithm.
560 
561    Input Parameters:
562    + tao - the Tao context
563    . func - function pointer for the regularizer value and gradient evaluation
564    - ctx - user context for the regularizer
565 
566    Level: advanced
567 @*/
568 PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (*func)(Tao,Vec,PetscReal *,Vec,void*),void *ctx)
569 {
570   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
571 
572   PetscFunctionBegin;
573   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
574   if (ctx) {
575     gn->reg_obj_ctx = ctx;
576   }
577   if (func) {
578     gn->regularizerobjandgrad = func;
579   }
580   PetscFunctionReturn(0);
581 }
582 
583 /*@C
584    TaoBRGNSetRegularizerHessianRoutine - Sets the user-defined regularizer call-back
585    function into the algorithm.
586 
587    Input Parameters:
588    + tao - the Tao context
589    . Hreg - user-created matrix for the Hessian of the regularization term
590    . func - function pointer for the regularizer Hessian evaluation
591    - ctx - user context for the regularizer Hessian
592 
593    Level: advanced
594 @*/
595 PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao,Mat Hreg,PetscErrorCode (*func)(Tao,Vec,Mat,void*),void *ctx)
596 {
597   TAO_BRGN       *gn = (TAO_BRGN *)tao->data;
598   PetscErrorCode ierr;
599 
600   PetscFunctionBegin;
601   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
602   if (Hreg) {
603     PetscValidHeaderSpecific(Hreg,MAT_CLASSID,2);
604     PetscCheckSameComm(tao,1,Hreg,2);
605   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONG,"NULL Hessian detected! User must provide valid Hessian for the regularizer.");
606   if (ctx) {
607     gn->reg_hess_ctx = ctx;
608   }
609   if (func) {
610     gn->regularizerhessian = func;
611   }
612   if (Hreg) {
613     ierr = PetscObjectReference((PetscObject)Hreg);CHKERRQ(ierr);
614     ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr);
615     gn->Hreg = Hreg;
616   }
617   PetscFunctionReturn(0);
618 }
619