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