xref: /petsc/src/tao/leastsquares/impls/brgn/brgn.c (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
1463fc0ecSAlp Dener #include <../src/tao/leastsquares/impls/brgn/brgn.h> /*I "petsctao.h" I*/
2737f463aSAlp Dener 
3470ec3f8SXiang Huang #define BRGN_REGULARIZATION_USER   0
4470ec3f8SXiang Huang #define BRGN_REGULARIZATION_L2PROX 1
5a1c74439SHansol Suh #define BRGN_REGULARIZATION_L2PURE 2
6a1c74439SHansol Suh #define BRGN_REGULARIZATION_L1DICT 3
7cd1c4666STristan Konolige #define BRGN_REGULARIZATION_LM     4
8cd1c4666STristan Konolige #define BRGN_REGULARIZATION_TYPES  5
9a3c390cfSAlp Dener 
10cd1c4666STristan Konolige static const char *BRGN_REGULARIZATION_TABLE[64] = {"user", "l2prox", "l2pure", "l1dict", "lm"};
11a3c390cfSAlp Dener 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode GNHessianProd(Mat H, Vec in, Vec out)
13d71ae5a4SJacob Faibussowitsch {
140d71dc2bSXiang Huang   TAO_BRGN *gn;
150d71dc2bSXiang Huang 
160d71dc2bSXiang Huang   PetscFunctionBegin;
179566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(H, &gn));
189566063dSJacob Faibussowitsch   PetscCall(MatMult(gn->subsolver->ls_jac, in, gn->r_work));
199566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(gn->subsolver->ls_jac, gn->r_work, out));
20a3c390cfSAlp Dener   switch (gn->reg_type) {
21470ec3f8SXiang Huang   case BRGN_REGULARIZATION_USER:
229566063dSJacob Faibussowitsch     PetscCall(MatMult(gn->Hreg, in, gn->x_work));
239566063dSJacob Faibussowitsch     PetscCall(VecAXPY(out, gn->lambda, gn->x_work));
24a3c390cfSAlp Dener     break;
25d71ae5a4SJacob Faibussowitsch   case BRGN_REGULARIZATION_L2PURE:
26d71ae5a4SJacob Faibussowitsch     PetscCall(VecAXPY(out, gn->lambda, in));
27d71ae5a4SJacob Faibussowitsch     break;
28d71ae5a4SJacob Faibussowitsch   case BRGN_REGULARIZATION_L2PROX:
29d71ae5a4SJacob Faibussowitsch     PetscCall(VecAXPY(out, gn->lambda, in));
30d71ae5a4SJacob Faibussowitsch     break;
31470ec3f8SXiang Huang   case BRGN_REGULARIZATION_L1DICT:
32a3c390cfSAlp Dener     /* out = out + lambda*D'*(diag.*(D*in)) */
33a3c390cfSAlp Dener     if (gn->D) {
349566063dSJacob Faibussowitsch       PetscCall(MatMult(gn->D, in, gn->y)); /* y = D*in */
35a3c390cfSAlp Dener     } else {
369566063dSJacob Faibussowitsch       PetscCall(VecCopy(in, gn->y));
37a3c390cfSAlp Dener     }
389566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(gn->y_work, gn->diag, gn->y)); /* y_work = diag.*(D*in), where diag = epsilon^2 ./ sqrt(x.^2+epsilon^2).^3 */
39a3c390cfSAlp Dener     if (gn->D) {
409566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(gn->D, gn->y_work, gn->x_work)); /* x_work = D'*(diag.*(D*in)) */
41a3c390cfSAlp Dener     } else {
429566063dSJacob Faibussowitsch       PetscCall(VecCopy(gn->y_work, gn->x_work));
43a3c390cfSAlp Dener     }
449566063dSJacob Faibussowitsch     PetscCall(VecAXPY(out, gn->lambda, gn->x_work));
45a3c390cfSAlp Dener     break;
46cd1c4666STristan Konolige   case BRGN_REGULARIZATION_LM:
479566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(gn->x_work, gn->damping, in));
489566063dSJacob Faibussowitsch     PetscCall(VecAXPY(out, 1, gn->x_work));
49cd1c4666STristan Konolige     break;
50a3c390cfSAlp Dener   }
513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
520d71dc2bSXiang Huang }
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeDamping(TAO_BRGN *gn)
54d71ae5a4SJacob Faibussowitsch {
55cd1c4666STristan Konolige   const PetscScalar *diag_ary;
56cd1c4666STristan Konolige   PetscScalar       *damping_ary;
57cd1c4666STristan Konolige   PetscInt           i, n;
58cd1c4666STristan Konolige 
59cd1c4666STristan Konolige   PetscFunctionBegin;
60cd1c4666STristan Konolige   /* update damping */
619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gn->damping, &damping_ary));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(gn->diag, &diag_ary));
639566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gn->damping, &n));
64ad540459SPierre Jolivet   for (i = 0; i < n; i++) damping_ary[i] = PetscClipInterval(diag_ary[i], PETSC_SQRT_MACHINE_EPSILON, PetscSqrtReal(PETSC_MAX_REAL));
659566063dSJacob Faibussowitsch   PetscCall(VecScale(gn->damping, gn->lambda));
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gn->damping, &damping_ary));
679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(gn->diag, &diag_ary));
683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69cd1c4666STristan Konolige }
70cd1c4666STristan Konolige 
71d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNGetDampingVector(Tao tao, Vec *d)
72d71ae5a4SJacob Faibussowitsch {
73cd1c4666STristan Konolige   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
74cd1c4666STristan Konolige 
75cd1c4666STristan Konolige   PetscFunctionBegin;
763c859ba3SBarry Smith   PetscCheck(gn->reg_type == BRGN_REGULARIZATION_LM, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Damping vector is only available if regularization type is lm.");
77cd1c4666STristan Konolige   *d = gn->damping;
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79cd1c4666STristan Konolige }
800d71dc2bSXiang Huang 
81d71ae5a4SJacob Faibussowitsch static PetscErrorCode GNObjectiveGradientEval(Tao tao, Vec X, PetscReal *fcn, Vec G, void *ptr)
82d71ae5a4SJacob Faibussowitsch {
830d71dc2bSXiang Huang   TAO_BRGN   *gn = (TAO_BRGN *)ptr;
848e85b1b3SXiang Huang   PetscInt    K; /* dimension of D*X */
857cea06e1SXiang Huang   PetscScalar yESum;
86a3c390cfSAlp Dener   PetscReal   f_reg;
870d71dc2bSXiang Huang 
880d71dc2bSXiang Huang   PetscFunctionBegin;
898e85b1b3SXiang Huang   /* compute objective *fcn*/
90a3c390cfSAlp Dener   /* compute first term 0.5*||ls_res||_2^2 */
919566063dSJacob Faibussowitsch   PetscCall(TaoComputeResidual(tao, X, tao->ls_res));
929566063dSJacob Faibussowitsch   PetscCall(VecDot(tao->ls_res, tao->ls_res, fcn));
93a3c390cfSAlp Dener   *fcn *= 0.5;
94a3c390cfSAlp Dener   /* compute gradient G */
959566063dSJacob Faibussowitsch   PetscCall(TaoComputeResidualJacobian(tao, X, tao->ls_jac, tao->ls_jac_pre));
969566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->ls_jac, tao->ls_res, G));
97a3c390cfSAlp Dener   /* add the regularization contribution */
98a3c390cfSAlp Dener   switch (gn->reg_type) {
99470ec3f8SXiang Huang   case BRGN_REGULARIZATION_USER:
1009566063dSJacob Faibussowitsch     PetscCall((*gn->regularizerobjandgrad)(tao, X, &f_reg, gn->x_work, gn->reg_obj_ctx));
101a3c390cfSAlp Dener     *fcn += gn->lambda * f_reg;
1029566063dSJacob Faibussowitsch     PetscCall(VecAXPY(G, gn->lambda, gn->x_work));
103a3c390cfSAlp Dener     break;
104a1c74439SHansol Suh   case BRGN_REGULARIZATION_L2PURE:
105a1c74439SHansol Suh     /* compute f = f + lambda*0.5*xk'*xk */
1069566063dSJacob Faibussowitsch     PetscCall(VecDot(X, X, &f_reg));
107a1c74439SHansol Suh     *fcn += gn->lambda * 0.5 * f_reg;
108a1c74439SHansol Suh     /* compute G = G + lambda*xk */
1099566063dSJacob Faibussowitsch     PetscCall(VecAXPY(G, gn->lambda, X));
110a1c74439SHansol Suh     break;
111470ec3f8SXiang Huang   case BRGN_REGULARIZATION_L2PROX:
1121fc140a9SXiang Huang     /* compute f = f + lambda*0.5*(xk - xkm1)'*(xk - xkm1) */
1139566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(gn->x_work, 1.0, -1.0, 0.0, X, gn->x_old));
1149566063dSJacob Faibussowitsch     PetscCall(VecDot(gn->x_work, gn->x_work, &f_reg));
115a3c390cfSAlp Dener     *fcn += gn->lambda * 0.5 * f_reg;
116a3c390cfSAlp Dener     /* compute G = G + lambda*(xk - xkm1) */
1179566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(G, gn->lambda, -gn->lambda, 1.0, X, gn->x_old));
118a3c390cfSAlp Dener     break;
119470ec3f8SXiang Huang   case BRGN_REGULARIZATION_L1DICT:
120a3c390cfSAlp Dener     /* compute f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x*/
121a3c390cfSAlp Dener     if (gn->D) {
1229566063dSJacob Faibussowitsch       PetscCall(MatMult(gn->D, X, gn->y)); /* y = D*x */
123a3c390cfSAlp Dener     } else {
1249566063dSJacob Faibussowitsch       PetscCall(VecCopy(X, gn->y));
125a3c390cfSAlp Dener     }
1269566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(gn->y_work, gn->y, gn->y));
1279566063dSJacob Faibussowitsch     PetscCall(VecShift(gn->y_work, gn->epsilon * gn->epsilon));
1289566063dSJacob Faibussowitsch     PetscCall(VecSqrtAbs(gn->y_work)); /* gn->y_work = sqrt(y.^2+epsilon^2) */
1299566063dSJacob Faibussowitsch     PetscCall(VecSum(gn->y_work, &yESum));
1309566063dSJacob Faibussowitsch     PetscCall(VecGetSize(gn->y, &K));
131a3c390cfSAlp Dener     *fcn += gn->lambda * (yESum - K * gn->epsilon);
1327cea06e1SXiang Huang     /* compute G = G + lambda*D'*(y./sqrt(y.^2+epsilon^2)),where y = D*x */
1339566063dSJacob Faibussowitsch     PetscCall(VecPointwiseDivide(gn->y_work, gn->y, gn->y_work)); /* reuse y_work = y./sqrt(y.^2+epsilon^2) */
134a3c390cfSAlp Dener     if (gn->D) {
1359566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(gn->D, gn->y_work, gn->x_work));
136a3c390cfSAlp Dener     } else {
1379566063dSJacob Faibussowitsch       PetscCall(VecCopy(gn->y_work, gn->x_work));
138a3c390cfSAlp Dener     }
1399566063dSJacob Faibussowitsch     PetscCall(VecAXPY(G, gn->lambda, gn->x_work));
140a3c390cfSAlp Dener     break;
141a3c390cfSAlp Dener   }
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1430d71dc2bSXiang Huang }
1440d71dc2bSXiang Huang 
145d71ae5a4SJacob Faibussowitsch static PetscErrorCode GNComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
146d71ae5a4SJacob Faibussowitsch {
1478ac80d48SXiang Huang   TAO_BRGN    *gn = (TAO_BRGN *)ptr;
148cd1c4666STristan Konolige   PetscInt     i, n, cstart, cend;
149cd1c4666STristan Konolige   PetscScalar *cnorms, *diag_ary;
150737f463aSAlp Dener 
151737f463aSAlp Dener   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   PetscCall(TaoComputeResidualJacobian(tao, X, tao->ls_jac, tao->ls_jac_pre));
15348a46eb9SPierre Jolivet   if (gn->mat_explicit) PetscCall(MatTransposeMatMult(tao->ls_jac, tao->ls_jac, MAT_REUSE_MATRIX, PETSC_DEFAULT, &gn->H));
1540d71dc2bSXiang Huang 
155a3c390cfSAlp Dener   switch (gn->reg_type) {
156470ec3f8SXiang Huang   case BRGN_REGULARIZATION_USER:
1579566063dSJacob Faibussowitsch     PetscCall((*gn->regularizerhessian)(tao, X, gn->Hreg, gn->reg_hess_ctx));
1581baa6e33SBarry Smith     if (gn->mat_explicit) PetscCall(MatAXPY(gn->H, 1.0, gn->Hreg, DIFFERENT_NONZERO_PATTERN));
159a3c390cfSAlp Dener     break;
160a1c74439SHansol Suh   case BRGN_REGULARIZATION_L2PURE:
1611baa6e33SBarry Smith     if (gn->mat_explicit) PetscCall(MatShift(gn->H, gn->lambda));
162a1c74439SHansol Suh     break;
163470ec3f8SXiang Huang   case BRGN_REGULARIZATION_L2PROX:
1641baa6e33SBarry Smith     if (gn->mat_explicit) PetscCall(MatShift(gn->H, gn->lambda));
165a3c390cfSAlp Dener     break;
166470ec3f8SXiang Huang   case BRGN_REGULARIZATION_L1DICT:
1677cea06e1SXiang Huang     /* 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 */
168a3c390cfSAlp Dener     if (gn->D) {
1699566063dSJacob Faibussowitsch       PetscCall(MatMult(gn->D, X, gn->y)); /* y = D*x */
170a3c390cfSAlp Dener     } else {
1719566063dSJacob Faibussowitsch       PetscCall(VecCopy(X, gn->y));
172a3c390cfSAlp Dener     }
1739566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(gn->y_work, gn->y, gn->y));
1749566063dSJacob Faibussowitsch     PetscCall(VecShift(gn->y_work, gn->epsilon * gn->epsilon));
1759566063dSJacob Faibussowitsch     PetscCall(VecCopy(gn->y_work, gn->diag));                    /* gn->diag = y.^2+epsilon^2 */
1769566063dSJacob Faibussowitsch     PetscCall(VecSqrtAbs(gn->y_work));                           /* gn->y_work = sqrt(y.^2+epsilon^2) */
1779566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(gn->diag, gn->y_work, gn->diag)); /* gn->diag = sqrt(y.^2+epsilon^2).^3 */
1789566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(gn->diag));
1799566063dSJacob Faibussowitsch     PetscCall(VecScale(gn->diag, gn->epsilon * gn->epsilon));
1801baa6e33SBarry Smith     if (gn->mat_explicit) PetscCall(MatDiagonalSet(gn->H, gn->diag, ADD_VALUES));
181a3c390cfSAlp Dener     break;
182cd1c4666STristan Konolige   case BRGN_REGULARIZATION_LM:
183cd1c4666STristan Konolige     /* compute diagonal of J^T J */
1849566063dSJacob Faibussowitsch     PetscCall(MatGetSize(gn->parent->ls_jac, NULL, &n));
1859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &cnorms));
1869566063dSJacob Faibussowitsch     PetscCall(MatGetColumnNorms(gn->parent->ls_jac, NORM_2, cnorms));
1879566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangeColumn(gn->parent->ls_jac, &cstart, &cend));
1889566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gn->diag, &diag_ary));
189ad540459SPierre Jolivet     for (i = 0; i < cend - cstart; i++) diag_ary[i] = cnorms[cstart + i] * cnorms[cstart + i];
1909566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gn->diag, &diag_ary));
1919566063dSJacob Faibussowitsch     PetscCall(PetscFree(cnorms));
1929566063dSJacob Faibussowitsch     PetscCall(ComputeDamping(gn));
1931baa6e33SBarry Smith     if (gn->mat_explicit) PetscCall(MatDiagonalSet(gn->H, gn->damping, ADD_VALUES));
194cd1c4666STristan Konolige     break;
195a3c390cfSAlp Dener   }
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197e1e80dc8SAlp Dener }
198e1e80dc8SAlp Dener 
199d71ae5a4SJacob Faibussowitsch static PetscErrorCode GNHookFunction(Tao tao, PetscInt iter, void *ctx)
200d71ae5a4SJacob Faibussowitsch {
2018fcddce6SStefano Zampini   TAO_BRGN *gn = (TAO_BRGN *)ctx;
202e1e80dc8SAlp Dener 
203e1e80dc8SAlp Dener   PetscFunctionBegin;
204e1e80dc8SAlp Dener   /* Update basic tao information from the subsolver */
205e1e80dc8SAlp Dener   gn->parent->nfuncs      = tao->nfuncs;
206e1e80dc8SAlp Dener   gn->parent->ngrads      = tao->ngrads;
207e1e80dc8SAlp Dener   gn->parent->nfuncgrads  = tao->nfuncgrads;
208e1e80dc8SAlp Dener   gn->parent->nhess       = tao->nhess;
209e1e80dc8SAlp Dener   gn->parent->niter       = tao->niter;
210e1e80dc8SAlp Dener   gn->parent->ksp_its     = tao->ksp_its;
211e1e80dc8SAlp Dener   gn->parent->ksp_tot_its = tao->ksp_tot_its;
212cd1c4666STristan Konolige   gn->parent->fc          = tao->fc;
2139566063dSJacob Faibussowitsch   PetscCall(TaoGetConvergedReason(tao, &gn->parent->reason));
214e1e80dc8SAlp Dener   /* Update the solution vectors */
215e1e80dc8SAlp Dener   if (iter == 0) {
2169566063dSJacob Faibussowitsch     PetscCall(VecSet(gn->x_old, 0.0));
217e1e80dc8SAlp Dener   } else {
2189566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, gn->x_old));
2199566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, gn->parent->solution));
220e1e80dc8SAlp Dener   }
221e1e80dc8SAlp Dener   /* Update the gradient */
2229566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, gn->parent->gradient));
223cd1c4666STristan Konolige 
224cd1c4666STristan Konolige   /* Update damping parameter for LM */
225cd1c4666STristan Konolige   if (gn->reg_type == BRGN_REGULARIZATION_LM) {
226cd1c4666STristan Konolige     if (iter > 0) {
227cd1c4666STristan Konolige       if (gn->fc_old > tao->fc) {
228cd1c4666STristan Konolige         gn->lambda = gn->lambda * gn->downhill_lambda_change;
229cd1c4666STristan Konolige       } else {
230cd1c4666STristan Konolige         /* uphill step */
231cd1c4666STristan Konolige         gn->lambda = gn->lambda * gn->uphill_lambda_change;
232cd1c4666STristan Konolige       }
233cd1c4666STristan Konolige     }
234cd1c4666STristan Konolige     gn->fc_old = tao->fc;
235cd1c4666STristan Konolige   }
236cd1c4666STristan Konolige 
237e1e80dc8SAlp Dener   /* Call general purpose update function */
2381baa6e33SBarry Smith   if (gn->parent->ops->update) PetscCall((*gn->parent->ops->update)(gn->parent, gn->parent->niter, gn->parent->user_update));
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240737f463aSAlp Dener }
241737f463aSAlp Dener 
242d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BRGN(Tao tao)
243d71ae5a4SJacob Faibussowitsch {
244737f463aSAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
245737f463aSAlp Dener 
246737f463aSAlp Dener   PetscFunctionBegin;
2479566063dSJacob Faibussowitsch   PetscCall(TaoSolve(gn->subsolver));
248e1e80dc8SAlp Dener   /* Update basic tao information from the subsolver */
249e1e80dc8SAlp Dener   tao->nfuncs      = gn->subsolver->nfuncs;
250e1e80dc8SAlp Dener   tao->ngrads      = gn->subsolver->ngrads;
251e1e80dc8SAlp Dener   tao->nfuncgrads  = gn->subsolver->nfuncgrads;
252e1e80dc8SAlp Dener   tao->nhess       = gn->subsolver->nhess;
253e1e80dc8SAlp Dener   tao->niter       = gn->subsolver->niter;
254e1e80dc8SAlp Dener   tao->ksp_its     = gn->subsolver->ksp_its;
255e1e80dc8SAlp Dener   tao->ksp_tot_its = gn->subsolver->ksp_tot_its;
2569566063dSJacob Faibussowitsch   PetscCall(TaoGetConvergedReason(gn->subsolver, &tao->reason));
257e1e80dc8SAlp Dener   /* Update vectors */
2589566063dSJacob Faibussowitsch   PetscCall(VecCopy(gn->subsolver->solution, tao->solution));
2599566063dSJacob Faibussowitsch   PetscCall(VecCopy(gn->subsolver->gradient, tao->gradient));
2603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
261737f463aSAlp Dener }
262737f463aSAlp Dener 
263d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_BRGN(Tao tao, PetscOptionItems *PetscOptionsObject)
264d71ae5a4SJacob Faibussowitsch {
265737f463aSAlp Dener   TAO_BRGN     *gn = (TAO_BRGN *)tao->data;
266cd1c4666STristan Konolige   TaoLineSearch ls;
267737f463aSAlp Dener 
268737f463aSAlp Dener   PetscFunctionBegin;
269d0609cedSBarry Smith   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.");
2709566063dSJacob Faibussowitsch   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));
2719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_brgn_regularizer_weight", "regularizer weight (default 1e-4)", "", gn->lambda, &gn->lambda, NULL));
2729566063dSJacob Faibussowitsch   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));
2739566063dSJacob Faibussowitsch   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));
2749566063dSJacob Faibussowitsch   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));
2759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_brgn_regularization_type", "regularization type", "", BRGN_REGULARIZATION_TABLE, BRGN_REGULARIZATION_TYPES, BRGN_REGULARIZATION_TABLE[gn->reg_type], &gn->reg_type, NULL));
276d0609cedSBarry Smith   PetscOptionsHeadEnd();
277cd1c4666STristan Konolige   /* set unit line search direction as the default when using the lm regularizer */
278cd1c4666STristan Konolige   if (gn->reg_type == BRGN_REGULARIZATION_LM) {
2799566063dSJacob Faibussowitsch     PetscCall(TaoGetLineSearch(gn->subsolver, &ls));
2809566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(ls, TAOLINESEARCHUNIT));
281cd1c4666STristan Konolige   }
2829566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(gn->subsolver));
2833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
284737f463aSAlp Dener }
285737f463aSAlp Dener 
286d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BRGN(Tao tao, PetscViewer viewer)
287d71ae5a4SJacob Faibussowitsch {
288737f463aSAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
289737f463aSAlp Dener 
290737f463aSAlp Dener   PetscFunctionBegin;
2919566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
2929566063dSJacob Faibussowitsch   PetscCall(TaoView(gn->subsolver, viewer));
2939566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
2943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
295737f463aSAlp Dener }
296737f463aSAlp Dener 
297d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BRGN(Tao tao)
298d71ae5a4SJacob Faibussowitsch {
299737f463aSAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
300737f463aSAlp Dener   PetscBool is_bnls, is_bntr, is_bntl;
3018e85b1b3SXiang Huang   PetscInt  i, n, N, K; /* dict has size K*N*/
302737f463aSAlp Dener 
303737f463aSAlp Dener   PetscFunctionBegin;
3043c859ba3SBarry Smith   PetscCheck(tao->ls_res, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetResidualRoutine() must be called before setup!");
3059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNLS, &is_bnls));
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNTR, &is_bntr));
3079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNTL, &is_bntl));
3083c859ba3SBarry Smith   PetscCheck((!is_bnls && !is_bntr && !is_bntl) || tao->ls_jac, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetResidualJacobianRoutine() must be called before setup!");
30948a46eb9SPierre Jolivet   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
31048a46eb9SPierre Jolivet   if (!gn->x_work) PetscCall(VecDuplicate(tao->solution, &gn->x_work));
31148a46eb9SPierre Jolivet   if (!gn->r_work) PetscCall(VecDuplicate(tao->ls_res, &gn->r_work));
312e1e80dc8SAlp Dener   if (!gn->x_old) {
3139566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &gn->x_old));
3149566063dSJacob Faibussowitsch     PetscCall(VecSet(gn->x_old, 0.0));
315e1e80dc8SAlp Dener   }
3167cea06e1SXiang Huang 
317470ec3f8SXiang Huang   if (BRGN_REGULARIZATION_L1DICT == gn->reg_type) {
3182036730cSSajid Ali     if (!gn->y) {
31930eeff36SXiang Huang       if (gn->D) {
3209566063dSJacob Faibussowitsch         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 */
3219566063dSJacob Faibussowitsch         PetscCall(MatCreateVecs(gn->D, NULL, &gn->y));
32230eeff36SXiang Huang       } else {
323*aaa8cc7dSPierre Jolivet         PetscCall(VecDuplicate(tao->solution, &gn->y)); /* If user does not setup dict matrix, use identity matrix, K=N */
32430eeff36SXiang Huang       }
3259566063dSJacob Faibussowitsch       PetscCall(VecSet(gn->y, 0.0));
3267cea06e1SXiang Huang     }
32748a46eb9SPierre Jolivet     if (!gn->y_work) PetscCall(VecDuplicate(gn->y, &gn->y_work));
3288ac80d48SXiang Huang     if (!gn->diag) {
3299566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(gn->y, &gn->diag));
3309566063dSJacob Faibussowitsch       PetscCall(VecSet(gn->diag, 0.0));
3318ac80d48SXiang Huang     }
33230eeff36SXiang Huang   }
333cd1c4666STristan Konolige   if (BRGN_REGULARIZATION_LM == gn->reg_type) {
33448a46eb9SPierre Jolivet     if (!gn->diag) PetscCall(MatCreateVecs(tao->ls_jac, &gn->diag, NULL));
33548a46eb9SPierre Jolivet     if (!gn->damping) PetscCall(MatCreateVecs(tao->ls_jac, &gn->damping, NULL));
336cd1c4666STristan Konolige   }
3370d71dc2bSXiang Huang 
338e1e80dc8SAlp Dener   if (!tao->setupcalled) {
339737f463aSAlp Dener     /* Hessian setup */
3405eb5f4d6SAlp Dener     if (gn->mat_explicit) {
3419566063dSJacob Faibussowitsch       PetscCall(TaoComputeResidualJacobian(tao, tao->solution, tao->ls_jac, tao->ls_jac_pre));
3429566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(tao->ls_jac, tao->ls_jac, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &gn->H));
3435eb5f4d6SAlp Dener     } else {
3449566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(tao->solution, &n));
3459566063dSJacob Faibussowitsch       PetscCall(VecGetSize(tao->solution, &N));
3469566063dSJacob Faibussowitsch       PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &gn->H));
3479566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(gn->H, n, n, N, N));
3489566063dSJacob Faibussowitsch       PetscCall(MatSetType(gn->H, MATSHELL));
3499566063dSJacob Faibussowitsch       PetscCall(MatSetOption(gn->H, MAT_SYMMETRIC, PETSC_TRUE));
3509566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(gn->H, MATOP_MULT, (void (*)(void))GNHessianProd));
3519566063dSJacob Faibussowitsch       PetscCall(MatShellSetContext(gn->H, gn));
3525eb5f4d6SAlp Dener     }
3539566063dSJacob Faibussowitsch     PetscCall(MatSetUp(gn->H));
354a5b23f4aSJose E. Roman     /* Subsolver setup,include initial vector and dictionary D */
3559566063dSJacob Faibussowitsch     PetscCall(TaoSetUpdate(gn->subsolver, GNHookFunction, gn));
3569566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(gn->subsolver, tao->solution));
3571baa6e33SBarry Smith     if (tao->bounded) PetscCall(TaoSetVariableBounds(gn->subsolver, tao->XL, tao->XU));
3589566063dSJacob Faibussowitsch     PetscCall(TaoSetResidualRoutine(gn->subsolver, tao->ls_res, tao->ops->computeresidual, tao->user_lsresP));
3599566063dSJacob Faibussowitsch     PetscCall(TaoSetJacobianResidualRoutine(gn->subsolver, tao->ls_jac, tao->ls_jac, tao->ops->computeresidualjacobian, tao->user_lsjacP));
3609566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(gn->subsolver, NULL, GNObjectiveGradientEval, gn));
3619566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(gn->subsolver, gn->H, gn->H, GNComputeHessian, gn));
362e1e80dc8SAlp Dener     /* Propagate some options down */
3639566063dSJacob Faibussowitsch     PetscCall(TaoSetTolerances(gn->subsolver, tao->gatol, tao->grtol, tao->gttol));
3649566063dSJacob Faibussowitsch     PetscCall(TaoSetMaximumIterations(gn->subsolver, tao->max_it));
3659566063dSJacob Faibussowitsch     PetscCall(TaoSetMaximumFunctionEvaluations(gn->subsolver, tao->max_funcs));
366737f463aSAlp Dener     for (i = 0; i < tao->numbermonitors; ++i) {
3679566063dSJacob Faibussowitsch       PetscCall(TaoSetMonitor(gn->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
3689566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(tao->monitorcontext[i])));
369737f463aSAlp Dener     }
3709566063dSJacob Faibussowitsch     PetscCall(TaoSetUp(gn->subsolver));
371e1e80dc8SAlp Dener   }
3723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
373737f463aSAlp Dener }
374737f463aSAlp Dener 
375d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BRGN(Tao tao)
376d71ae5a4SJacob Faibussowitsch {
377737f463aSAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
378737f463aSAlp Dener 
379737f463aSAlp Dener   PetscFunctionBegin;
380737f463aSAlp Dener   if (tao->setupcalled) {
3819566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tao->gradient));
3829566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&gn->x_work));
3839566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&gn->r_work));
3849566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&gn->x_old));
3859566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&gn->diag));
3869566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&gn->y));
3879566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&gn->y_work));
388737f463aSAlp Dener   }
3899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gn->damping));
3909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gn->diag));
3919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&gn->H));
3929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&gn->D));
3939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&gn->Hreg));
3949566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&gn->subsolver));
395e1e80dc8SAlp Dener   gn->parent = NULL;
3969566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
398737f463aSAlp Dener }
399737f463aSAlp Dener 
4003850be85SAlp Dener /*MC
4013850be85SAlp Dener   TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares
4022fe279fdSBarry Smith             problems with bound constraints. This algorithm is a thin wrapper around `TAOBNTL`
403463fc0ecSAlp Dener             that constructs the Gauss-Newton problem with the user-provided least-squares
40460bb7533SHansol Suh             residual and Jacobian. The algorithm offers an L2-norm ("l2pure"), L2-norm proximal point ("l2prox")
40560bb7533SHansol Suh             regularizer, and L1-norm dictionary regularizer ("l1dict"), where we approximate the
40601b716f5SXiang Huang             L1-norm ||x||_1 by sum_i(sqrt(x_i^2+epsilon^2)-epsilon) with a small positive number epsilon.
407cd1c4666STristan Konolige             Also offered is the "lm" regularizer which uses a scaled diagonal of J^T J.
4082fe279fdSBarry Smith             With the "lm" regularizer, `TAOBRGN` is a Levenberg-Marquardt optimizer.
40901b716f5SXiang Huang             The user can also provide own regularization function.
4103850be85SAlp Dener 
4113850be85SAlp Dener   Options Database Keys:
412cd1c4666STristan Konolige + -tao_brgn_regularization_type - regularization type ("user", "l2prox", "l2pure", "l1dict", "lm") (default "l2prox")
413c061e8e2SXiang Huang . -tao_brgn_regularizer_weight  - regularizer weight (default 1e-4)
414c061e8e2SXiang Huang - -tao_brgn_l1_smooth_epsilon   - L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6)
4153850be85SAlp Dener 
4163850be85SAlp Dener   Level: beginner
4172fe279fdSBarry Smith 
4182fe279fdSBarry Smith .seealso: `Tao`, `TaoBRGNGetSubsolver()`, `TaoBRGNSetRegularizerWeight()`, `TaoBRGNSetL1SmoothEpsilon()`, `TaoBRGNSetDictionaryMatrix()`,
4192fe279fdSBarry Smith           `TaoBRGNSetRegularizerObjectiveAndGradientRoutine()`, `TaoBRGNSetRegularizerHessianRoutine()`
4203850be85SAlp Dener M*/
421d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao)
422d71ae5a4SJacob Faibussowitsch {
423737f463aSAlp Dener   TAO_BRGN *gn;
424737f463aSAlp Dener 
425737f463aSAlp Dener   PetscFunctionBegin;
4264dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&gn));
427737f463aSAlp Dener 
428737f463aSAlp Dener   tao->ops->destroy        = TaoDestroy_BRGN;
429737f463aSAlp Dener   tao->ops->setup          = TaoSetUp_BRGN;
430737f463aSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BRGN;
431737f463aSAlp Dener   tao->ops->view           = TaoView_BRGN;
432737f463aSAlp Dener   tao->ops->solve          = TaoSolve_BRGN;
433737f463aSAlp Dener 
4343ec1f749SStefano Zampini   tao->data                  = gn;
435d8bf7057SXiang Huang   gn->reg_type               = BRGN_REGULARIZATION_L2PROX;
436e1e80dc8SAlp Dener   gn->lambda                 = 1e-4;
4378ac80d48SXiang Huang   gn->epsilon                = 1e-6;
438cd1c4666STristan Konolige   gn->downhill_lambda_change = 1. / 5.;
439cd1c4666STristan Konolige   gn->uphill_lambda_change   = 1.5;
440e1e80dc8SAlp Dener   gn->parent                 = tao;
441737f463aSAlp Dener 
4429566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &gn->subsolver));
4439566063dSJacob Faibussowitsch   PetscCall(TaoSetType(gn->subsolver, TAOBNLS));
4449566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(gn->subsolver, "tao_brgn_subsolver_"));
4453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
446e1e80dc8SAlp Dener }
447e1e80dc8SAlp Dener 
448463fc0ecSAlp Dener /*@
4492fe279fdSBarry Smith   TaoBRGNGetSubsolver - Get the pointer to the subsolver inside a `TAOBRGN`
450e1e80dc8SAlp Dener 
451c3339decSBarry Smith   Collective
452e1e80dc8SAlp Dener 
453e1e80dc8SAlp Dener   Input Parameters:
454e1e80dc8SAlp Dener +  tao - the Tao solver context
4552fe279fdSBarry Smith -  subsolver - the `Tao` sub-solver context
4562fe279fdSBarry Smith 
4572fe279fdSBarry Smith   Level: advanced
4582fe279fdSBarry Smith 
4592fe279fdSBarry Smith .seealso: `Tao`, `Mat`, `TAOBRGN`
460e1e80dc8SAlp Dener @*/
461d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNGetSubsolver(Tao tao, Tao *subsolver)
462d71ae5a4SJacob Faibussowitsch {
463e1e80dc8SAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
464e1e80dc8SAlp Dener 
465e1e80dc8SAlp Dener   PetscFunctionBegin;
466e1e80dc8SAlp Dener   *subsolver = gn->subsolver;
4673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
468737f463aSAlp Dener }
469737f463aSAlp Dener 
470463fc0ecSAlp Dener /*@
471463fc0ecSAlp Dener   TaoBRGNSetRegularizerWeight - Set the regularizer weight for the Gauss-Newton least-squares algorithm
472737f463aSAlp Dener 
473c3339decSBarry Smith   Collective
474737f463aSAlp Dener 
475737f463aSAlp Dener   Input Parameters:
4762fe279fdSBarry Smith +  tao - the `Tao` solver context
4778e85b1b3SXiang Huang -  lambda - L1-norm regularizer weight
478463fc0ecSAlp Dener 
479463fc0ecSAlp Dener   Level: beginner
4802fe279fdSBarry Smith 
4812fe279fdSBarry Smith .seealso: `Tao`, `Mat`, `TAOBRGN`
482737f463aSAlp Dener @*/
483d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao, PetscReal lambda)
484d71ae5a4SJacob Faibussowitsch {
485737f463aSAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
486737f463aSAlp Dener 
4878ac80d48SXiang Huang   /* Initialize lambda here */
4880d71dc2bSXiang Huang 
489737f463aSAlp Dener   PetscFunctionBegin;
490737f463aSAlp Dener   gn->lambda = lambda;
4913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
492737f463aSAlp Dener }
4930d71dc2bSXiang Huang 
494463fc0ecSAlp Dener /*@
4958ac80d48SXiang Huang   TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm
4968ac80d48SXiang Huang 
497c3339decSBarry Smith   Collective
4988ac80d48SXiang Huang 
4998ac80d48SXiang Huang   Input Parameters:
5002fe279fdSBarry Smith +  tao - the `Tao` solver context
5018ac80d48SXiang Huang -  epsilon - L1-norm smooth approximation parameter
502463fc0ecSAlp Dener 
503463fc0ecSAlp Dener   Level: advanced
5042fe279fdSBarry Smith 
5052fe279fdSBarry Smith .seealso: `Tao`, `Mat`, `TAOBRGN`
5068ac80d48SXiang Huang @*/
507d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao, PetscReal epsilon)
508d71ae5a4SJacob Faibussowitsch {
5098ac80d48SXiang Huang   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
5108ac80d48SXiang Huang 
5118ac80d48SXiang Huang   /* Initialize epsilon here */
5128ac80d48SXiang Huang 
5138ac80d48SXiang Huang   PetscFunctionBegin;
5148ac80d48SXiang Huang   gn->epsilon = epsilon;
5153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5168ac80d48SXiang Huang }
5178e85b1b3SXiang Huang 
518463fc0ecSAlp Dener /*@
5198e85b1b3SXiang Huang    TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem)
5208e85b1b3SXiang Huang 
5218e85b1b3SXiang Huang    Input Parameters:
5222fe279fdSBarry Smith +  tao  - the `Tao` context
5232fe279fdSBarry Smith -  dict - the user specified dictionary matrix.  We allow to set a `NULL` dictionary, which means identity matrix by default
5248e85b1b3SXiang Huang 
525463fc0ecSAlp Dener     Level: advanced
5262fe279fdSBarry Smith 
5272fe279fdSBarry Smith .seealso: `Tao`, `Mat`, `TAOBRGN`
5288e85b1b3SXiang Huang @*/
529d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao, Mat dict)
530d71ae5a4SJacob Faibussowitsch {
5318e85b1b3SXiang Huang   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
5328e85b1b3SXiang Huang   PetscFunctionBegin;
5338e85b1b3SXiang Huang   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
5348e85b1b3SXiang Huang   if (dict) {
5358e85b1b3SXiang Huang     PetscValidHeaderSpecific(dict, MAT_CLASSID, 2);
536a3c390cfSAlp Dener     PetscCheckSameComm(tao, 1, dict, 2);
5379566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dict));
5388e85b1b3SXiang Huang   }
5399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&gn->D));
5401fc140a9SXiang Huang   gn->D = dict;
5413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5428e85b1b3SXiang Huang }
5438e85b1b3SXiang Huang 
544a3c390cfSAlp Dener /*@C
545463fc0ecSAlp Dener    TaoBRGNSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back
546463fc0ecSAlp Dener    function into the algorithm.
547463fc0ecSAlp Dener 
548463fc0ecSAlp Dener    Input Parameters:
549463fc0ecSAlp Dener + tao - the Tao context
550463fc0ecSAlp Dener . func - function pointer for the regularizer value and gradient evaluation
551463fc0ecSAlp Dener - ctx - user context for the regularizer
552463fc0ecSAlp Dener 
5532fe279fdSBarry Smith    Calling Sequence of `func`:
5542fe279fdSBarry Smith $  PetscErrorCode (*func)(Tao tao, Vec u, PetscReal val, Vec g, void *ctx)
5552fe279fdSBarry Smith +  tao - the `Tao` context
5562fe279fdSBarry Smith .  u - the location at which to compute the objective and gradient
5572fe279fdSBarry Smith .  val - location to store objective function value
5582fe279fdSBarry Smith .  g - location to store gradient
5592fe279fdSBarry Smith -  ctx - user context for the regularizer Hessian
5602fe279fdSBarry Smith 
561463fc0ecSAlp Dener    Level: advanced
5622fe279fdSBarry Smith 
5632fe279fdSBarry Smith .seealso: `Tao`, `Mat`, `TAOBRGN`
564a3c390cfSAlp Dener @*/
565d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), void *ctx)
566d71ae5a4SJacob Faibussowitsch {
567a3c390cfSAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
568a3c390cfSAlp Dener 
569a3c390cfSAlp Dener   PetscFunctionBegin;
570a3c390cfSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
571ad540459SPierre Jolivet   if (ctx) gn->reg_obj_ctx = ctx;
572ad540459SPierre Jolivet   if (func) gn->regularizerobjandgrad = func;
5733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
574a3c390cfSAlp Dener }
575a3c390cfSAlp Dener 
576a3c390cfSAlp Dener /*@C
577463fc0ecSAlp Dener    TaoBRGNSetRegularizerHessianRoutine - Sets the user-defined regularizer call-back
578463fc0ecSAlp Dener    function into the algorithm.
579463fc0ecSAlp Dener 
580463fc0ecSAlp Dener    Input Parameters:
5812fe279fdSBarry Smith +  tao - the `Tao` context
582463fc0ecSAlp Dener .  Hreg - user-created matrix for the Hessian of the regularization term
583463fc0ecSAlp Dener .  func - function pointer for the regularizer Hessian evaluation
584463fc0ecSAlp Dener -  ctx - user context for the regularizer Hessian
585463fc0ecSAlp Dener 
5862fe279fdSBarry Smith    Calling Sequence of `func`:
5872fe279fdSBarry Smith $  PetscErrorCode (*func)(Tao tao, Vec u, Mat H, void *ctx)
5882fe279fdSBarry Smith +  tao - the `Tao` context
5892fe279fdSBarry Smith .  u - the location at which to compute the Hessian
5902fe279fdSBarry Smith .  Hreg - user-created matrix for the Hessian of the regularization term
5912fe279fdSBarry Smith -  ctx - user context for the regularizer Hessian
5922fe279fdSBarry Smith 
593463fc0ecSAlp Dener    Level: advanced
5942fe279fdSBarry Smith 
5952fe279fdSBarry Smith   .seealso: `Tao`, `Mat`, `TAOBRGN`
596a3c390cfSAlp Dener @*/
597d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao, Mat Hreg, PetscErrorCode (*func)(Tao, Vec, Mat, void *), void *ctx)
598d71ae5a4SJacob Faibussowitsch {
599a3c390cfSAlp Dener   TAO_BRGN *gn = (TAO_BRGN *)tao->data;
600a3c390cfSAlp Dener 
601a3c390cfSAlp Dener   PetscFunctionBegin;
602a3c390cfSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
603a3c390cfSAlp Dener   if (Hreg) {
604a3c390cfSAlp Dener     PetscValidHeaderSpecific(Hreg, MAT_CLASSID, 2);
605a3c390cfSAlp Dener     PetscCheckSameComm(tao, 1, Hreg, 2);
60601b716f5SXiang Huang   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "NULL Hessian detected! User must provide valid Hessian for the regularizer.");
607ad540459SPierre Jolivet   if (ctx) gn->reg_hess_ctx = ctx;
608ad540459SPierre Jolivet   if (func) gn->regularizerhessian = func;
609a3c390cfSAlp Dener   if (Hreg) {
6109566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Hreg));
6119566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&gn->Hreg));
612a3c390cfSAlp Dener     gn->Hreg = Hreg;
613a3c390cfSAlp Dener   }
6143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
615a3c390cfSAlp Dener }
616