1463fc0ecSAlp Dener #include <../src/tao/leastsquares/impls/brgn/brgn.h> /*I "petsctao.h" I*/
2737f463aSAlp Dener
3c0b7dd19SHansol Suh const char *const TaoBRGNRegularizationTypes[] = {"user", "l2prox", "l2pure", "l1dict", "lm", "TaoBRGNRegularizationType", "TAOBRGN_REGULARIZATION_", NULL};
4a3c390cfSAlp Dener
GNHessianProd(Mat H,Vec in,Vec out)5d71ae5a4SJacob Faibussowitsch static PetscErrorCode GNHessianProd(Mat H, Vec in, Vec out)
6d71ae5a4SJacob Faibussowitsch {
70d71dc2bSXiang Huang TAO_BRGN *gn;
80d71dc2bSXiang Huang
90d71dc2bSXiang Huang PetscFunctionBegin;
109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(H, &gn));
119566063dSJacob Faibussowitsch PetscCall(MatMult(gn->subsolver->ls_jac, in, gn->r_work));
129566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(gn->subsolver->ls_jac, gn->r_work, out));
13a3c390cfSAlp Dener switch (gn->reg_type) {
14c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_USER:
159566063dSJacob Faibussowitsch PetscCall(MatMult(gn->Hreg, in, gn->x_work));
169566063dSJacob Faibussowitsch PetscCall(VecAXPY(out, gn->lambda, gn->x_work));
17a3c390cfSAlp Dener break;
18c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PURE:
19d71ae5a4SJacob Faibussowitsch PetscCall(VecAXPY(out, gn->lambda, in));
20d71ae5a4SJacob Faibussowitsch break;
21c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PROX:
22d71ae5a4SJacob Faibussowitsch PetscCall(VecAXPY(out, gn->lambda, in));
23d71ae5a4SJacob Faibussowitsch break;
24c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L1DICT:
25a3c390cfSAlp Dener /* out = out + lambda*D'*(diag.*(D*in)) */
26a3c390cfSAlp Dener if (gn->D) {
279566063dSJacob Faibussowitsch PetscCall(MatMult(gn->D, in, gn->y)); /* y = D*in */
28a3c390cfSAlp Dener } else {
299566063dSJacob Faibussowitsch PetscCall(VecCopy(in, gn->y));
30a3c390cfSAlp Dener }
319566063dSJacob 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 */
32a3c390cfSAlp Dener if (gn->D) {
339566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(gn->D, gn->y_work, gn->x_work)); /* x_work = D'*(diag.*(D*in)) */
34a3c390cfSAlp Dener } else {
359566063dSJacob Faibussowitsch PetscCall(VecCopy(gn->y_work, gn->x_work));
36a3c390cfSAlp Dener }
379566063dSJacob Faibussowitsch PetscCall(VecAXPY(out, gn->lambda, gn->x_work));
38a3c390cfSAlp Dener break;
39c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_LM:
409566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(gn->x_work, gn->damping, in));
419566063dSJacob Faibussowitsch PetscCall(VecAXPY(out, 1, gn->x_work));
42cd1c4666STristan Konolige break;
43a3c390cfSAlp Dener }
443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
450d71dc2bSXiang Huang }
ComputeDamping(TAO_BRGN * gn)46d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeDamping(TAO_BRGN *gn)
47d71ae5a4SJacob Faibussowitsch {
48cd1c4666STristan Konolige const PetscScalar *diag_ary;
49cd1c4666STristan Konolige PetscScalar *damping_ary;
50cd1c4666STristan Konolige PetscInt i, n;
51cd1c4666STristan Konolige
52cd1c4666STristan Konolige PetscFunctionBegin;
53cd1c4666STristan Konolige /* update damping */
549566063dSJacob Faibussowitsch PetscCall(VecGetArray(gn->damping, &damping_ary));
559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gn->diag, &diag_ary));
569566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gn->damping, &n));
57ad540459SPierre Jolivet for (i = 0; i < n; i++) damping_ary[i] = PetscClipInterval(diag_ary[i], PETSC_SQRT_MACHINE_EPSILON, PetscSqrtReal(PETSC_MAX_REAL));
589566063dSJacob Faibussowitsch PetscCall(VecScale(gn->damping, gn->lambda));
599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gn->damping, &damping_ary));
609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gn->diag, &diag_ary));
613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
62cd1c4666STristan Konolige }
63c0b7dd19SHansol Suh /*@
64c0b7dd19SHansol Suh TaoBRGNGetDampingVector - Get the damping vector $\mathrm{diag}(J^T J)$ from a `TAOBRGN` with `TAOBRGN_REGULARIZATION_LM` regularization
65cd1c4666STristan Konolige
66c0b7dd19SHansol Suh Collective
67c0b7dd19SHansol Suh
68c0b7dd19SHansol Suh Input Parameter:
69c0b7dd19SHansol Suh . tao - a `Tao` of type `TAOBRGN` with `TAOBRGN_REGULARIZATION_LM` regularization
70c0b7dd19SHansol Suh
71c0b7dd19SHansol Suh Output Parameter:
72c0b7dd19SHansol Suh . d - the damping vector
73c0b7dd19SHansol Suh
74c0b7dd19SHansol Suh Level: developer
75c0b7dd19SHansol Suh
76c0b7dd19SHansol Suh .seealso: [](ch_tao), `Tao`, `TAOBRGN`, `TaoBRGNRegularzationTypes`
77c0b7dd19SHansol Suh @*/
TaoBRGNGetDampingVector(Tao tao,Vec * d)78d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBRGNGetDampingVector(Tao tao, Vec *d)
79d71ae5a4SJacob Faibussowitsch {
80c0b7dd19SHansol Suh PetscFunctionBegin;
81c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
82c0b7dd19SHansol Suh PetscAssertPointer(d, 2);
83c0b7dd19SHansol Suh PetscUseMethod((PetscObject)tao, "TaoBRGNGetDampingVector_C", (Tao, Vec *), (tao, d));
84c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
85c0b7dd19SHansol Suh }
86c0b7dd19SHansol Suh
TaoBRGNGetDampingVector_BRGN(Tao tao,Vec * d)87c0b7dd19SHansol Suh static PetscErrorCode TaoBRGNGetDampingVector_BRGN(Tao tao, Vec *d)
88c0b7dd19SHansol Suh {
89cd1c4666STristan Konolige TAO_BRGN *gn = (TAO_BRGN *)tao->data;
90cd1c4666STristan Konolige
91cd1c4666STristan Konolige PetscFunctionBegin;
92c0b7dd19SHansol Suh PetscCheck(gn->reg_type == TAOBRGN_REGULARIZATION_LM, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Damping vector is only available if regularization type is lm.");
93cd1c4666STristan Konolige *d = gn->damping;
943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
95cd1c4666STristan Konolige }
960d71dc2bSXiang Huang
GNObjectiveGradientEval(Tao tao,Vec X,PetscReal * fcn,Vec G,void * ptr)97d71ae5a4SJacob Faibussowitsch static PetscErrorCode GNObjectiveGradientEval(Tao tao, Vec X, PetscReal *fcn, Vec G, void *ptr)
98d71ae5a4SJacob Faibussowitsch {
990d71dc2bSXiang Huang TAO_BRGN *gn = (TAO_BRGN *)ptr;
1008e85b1b3SXiang Huang PetscInt K; /* dimension of D*X */
1017cea06e1SXiang Huang PetscScalar yESum;
102a3c390cfSAlp Dener PetscReal f_reg;
1030d71dc2bSXiang Huang
1040d71dc2bSXiang Huang PetscFunctionBegin;
1058e85b1b3SXiang Huang /* compute objective *fcn*/
106a3c390cfSAlp Dener /* compute first term 0.5*||ls_res||_2^2 */
1079566063dSJacob Faibussowitsch PetscCall(TaoComputeResidual(tao, X, tao->ls_res));
1089566063dSJacob Faibussowitsch PetscCall(VecDot(tao->ls_res, tao->ls_res, fcn));
109a3c390cfSAlp Dener *fcn *= 0.5;
110a3c390cfSAlp Dener /* compute gradient G */
1119566063dSJacob Faibussowitsch PetscCall(TaoComputeResidualJacobian(tao, X, tao->ls_jac, tao->ls_jac_pre));
1129566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->ls_jac, tao->ls_res, G));
113a3c390cfSAlp Dener /* add the regularization contribution */
114a3c390cfSAlp Dener switch (gn->reg_type) {
115c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_USER:
1169566063dSJacob Faibussowitsch PetscCall((*gn->regularizerobjandgrad)(tao, X, &f_reg, gn->x_work, gn->reg_obj_ctx));
117a3c390cfSAlp Dener *fcn += gn->lambda * f_reg;
1189566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, gn->lambda, gn->x_work));
119a3c390cfSAlp Dener break;
120c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PURE:
121a1c74439SHansol Suh /* compute f = f + lambda*0.5*xk'*xk */
1229566063dSJacob Faibussowitsch PetscCall(VecDot(X, X, &f_reg));
123a1c74439SHansol Suh *fcn += gn->lambda * 0.5 * f_reg;
124a1c74439SHansol Suh /* compute G = G + lambda*xk */
1259566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, gn->lambda, X));
126a1c74439SHansol Suh break;
127c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PROX:
1281fc140a9SXiang Huang /* compute f = f + lambda*0.5*(xk - xkm1)'*(xk - xkm1) */
1299566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(gn->x_work, 1.0, -1.0, 0.0, X, gn->x_old));
1309566063dSJacob Faibussowitsch PetscCall(VecDot(gn->x_work, gn->x_work, &f_reg));
131a3c390cfSAlp Dener *fcn += gn->lambda * 0.5 * f_reg;
132a3c390cfSAlp Dener /* compute G = G + lambda*(xk - xkm1) */
1339566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(G, gn->lambda, -gn->lambda, 1.0, X, gn->x_old));
134a3c390cfSAlp Dener break;
135c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L1DICT:
136a3c390cfSAlp Dener /* compute f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x*/
137a3c390cfSAlp Dener if (gn->D) {
1389566063dSJacob Faibussowitsch PetscCall(MatMult(gn->D, X, gn->y)); /* y = D*x */
139a3c390cfSAlp Dener } else {
1409566063dSJacob Faibussowitsch PetscCall(VecCopy(X, gn->y));
141a3c390cfSAlp Dener }
1429566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(gn->y_work, gn->y, gn->y));
1439566063dSJacob Faibussowitsch PetscCall(VecShift(gn->y_work, gn->epsilon * gn->epsilon));
1449566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(gn->y_work)); /* gn->y_work = sqrt(y.^2+epsilon^2) */
1459566063dSJacob Faibussowitsch PetscCall(VecSum(gn->y_work, &yESum));
1469566063dSJacob Faibussowitsch PetscCall(VecGetSize(gn->y, &K));
147a3c390cfSAlp Dener *fcn += gn->lambda * (yESum - K * gn->epsilon);
1487cea06e1SXiang Huang /* compute G = G + lambda*D'*(y./sqrt(y.^2+epsilon^2)),where y = D*x */
1499566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(gn->y_work, gn->y, gn->y_work)); /* reuse y_work = y./sqrt(y.^2+epsilon^2) */
150a3c390cfSAlp Dener if (gn->D) {
1519566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(gn->D, gn->y_work, gn->x_work));
152a3c390cfSAlp Dener } else {
1539566063dSJacob Faibussowitsch PetscCall(VecCopy(gn->y_work, gn->x_work));
154a3c390cfSAlp Dener }
1559566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, gn->lambda, gn->x_work));
156a3c390cfSAlp Dener break;
157c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_LM:
158c0b7dd19SHansol Suh break;
159c0b7dd19SHansol Suh default:
160c0b7dd19SHansol Suh break;
161a3c390cfSAlp Dener }
1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1630d71dc2bSXiang Huang }
1640d71dc2bSXiang Huang
GNComputeHessian(Tao tao,Vec X,Mat H,Mat Hpre,void * ptr)165d71ae5a4SJacob Faibussowitsch static PetscErrorCode GNComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
166d71ae5a4SJacob Faibussowitsch {
1678ac80d48SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)ptr;
168cd1c4666STristan Konolige PetscInt i, n, cstart, cend;
169cd1c4666STristan Konolige PetscScalar *cnorms, *diag_ary;
170737f463aSAlp Dener
171737f463aSAlp Dener PetscFunctionBegin;
1729566063dSJacob Faibussowitsch PetscCall(TaoComputeResidualJacobian(tao, X, tao->ls_jac, tao->ls_jac_pre));
173fb842aefSJose E. Roman if (gn->mat_explicit) PetscCall(MatTransposeMatMult(tao->ls_jac, tao->ls_jac, MAT_REUSE_MATRIX, PETSC_DETERMINE, &gn->H));
1740d71dc2bSXiang Huang
175a3c390cfSAlp Dener switch (gn->reg_type) {
176c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_USER:
1779566063dSJacob Faibussowitsch PetscCall((*gn->regularizerhessian)(tao, X, gn->Hreg, gn->reg_hess_ctx));
1781baa6e33SBarry Smith if (gn->mat_explicit) PetscCall(MatAXPY(gn->H, 1.0, gn->Hreg, DIFFERENT_NONZERO_PATTERN));
179a3c390cfSAlp Dener break;
180c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PURE:
1811baa6e33SBarry Smith if (gn->mat_explicit) PetscCall(MatShift(gn->H, gn->lambda));
182a1c74439SHansol Suh break;
183c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PROX:
1841baa6e33SBarry Smith if (gn->mat_explicit) PetscCall(MatShift(gn->H, gn->lambda));
185a3c390cfSAlp Dener break;
186c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L1DICT:
1877cea06e1SXiang 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 */
188a3c390cfSAlp Dener if (gn->D) {
1899566063dSJacob Faibussowitsch PetscCall(MatMult(gn->D, X, gn->y)); /* y = D*x */
190a3c390cfSAlp Dener } else {
1919566063dSJacob Faibussowitsch PetscCall(VecCopy(X, gn->y));
192a3c390cfSAlp Dener }
1939566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(gn->y_work, gn->y, gn->y));
1949566063dSJacob Faibussowitsch PetscCall(VecShift(gn->y_work, gn->epsilon * gn->epsilon));
1959566063dSJacob Faibussowitsch PetscCall(VecCopy(gn->y_work, gn->diag)); /* gn->diag = y.^2+epsilon^2 */
1969566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(gn->y_work)); /* gn->y_work = sqrt(y.^2+epsilon^2) */
1979566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(gn->diag, gn->y_work, gn->diag)); /* gn->diag = sqrt(y.^2+epsilon^2).^3 */
1989566063dSJacob Faibussowitsch PetscCall(VecReciprocal(gn->diag));
1999566063dSJacob Faibussowitsch PetscCall(VecScale(gn->diag, gn->epsilon * gn->epsilon));
2001baa6e33SBarry Smith if (gn->mat_explicit) PetscCall(MatDiagonalSet(gn->H, gn->diag, ADD_VALUES));
201a3c390cfSAlp Dener break;
202c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_LM:
203cd1c4666STristan Konolige /* compute diagonal of J^T J */
2049566063dSJacob Faibussowitsch PetscCall(MatGetSize(gn->parent->ls_jac, NULL, &n));
2059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &cnorms));
2069566063dSJacob Faibussowitsch PetscCall(MatGetColumnNorms(gn->parent->ls_jac, NORM_2, cnorms));
2079566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(gn->parent->ls_jac, &cstart, &cend));
2089566063dSJacob Faibussowitsch PetscCall(VecGetArray(gn->diag, &diag_ary));
209ad540459SPierre Jolivet for (i = 0; i < cend - cstart; i++) diag_ary[i] = cnorms[cstart + i] * cnorms[cstart + i];
2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gn->diag, &diag_ary));
2119566063dSJacob Faibussowitsch PetscCall(PetscFree(cnorms));
2129566063dSJacob Faibussowitsch PetscCall(ComputeDamping(gn));
2131baa6e33SBarry Smith if (gn->mat_explicit) PetscCall(MatDiagonalSet(gn->H, gn->damping, ADD_VALUES));
214cd1c4666STristan Konolige break;
215c0b7dd19SHansol Suh default:
216c0b7dd19SHansol Suh break;
217a3c390cfSAlp Dener }
2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
219e1e80dc8SAlp Dener }
220e1e80dc8SAlp Dener
GNHookFunction(Tao tao,PetscInt iter,PetscCtx ctx)2212a8381b2SBarry Smith static PetscErrorCode GNHookFunction(Tao tao, PetscInt iter, PetscCtx ctx)
222d71ae5a4SJacob Faibussowitsch {
2238fcddce6SStefano Zampini TAO_BRGN *gn = (TAO_BRGN *)ctx;
224e1e80dc8SAlp Dener
225e1e80dc8SAlp Dener PetscFunctionBegin;
226e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */
227e1e80dc8SAlp Dener gn->parent->nfuncs = tao->nfuncs;
228e1e80dc8SAlp Dener gn->parent->ngrads = tao->ngrads;
229e1e80dc8SAlp Dener gn->parent->nfuncgrads = tao->nfuncgrads;
230e1e80dc8SAlp Dener gn->parent->nhess = tao->nhess;
231e1e80dc8SAlp Dener gn->parent->niter = tao->niter;
232e1e80dc8SAlp Dener gn->parent->ksp_its = tao->ksp_its;
233e1e80dc8SAlp Dener gn->parent->ksp_tot_its = tao->ksp_tot_its;
234cd1c4666STristan Konolige gn->parent->fc = tao->fc;
2359566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao, &gn->parent->reason));
236e1e80dc8SAlp Dener /* Update the solution vectors */
237e1e80dc8SAlp Dener if (iter == 0) {
2389566063dSJacob Faibussowitsch PetscCall(VecSet(gn->x_old, 0.0));
239e1e80dc8SAlp Dener } else {
2409566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, gn->x_old));
2419566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, gn->parent->solution));
242e1e80dc8SAlp Dener }
243e1e80dc8SAlp Dener /* Update the gradient */
2449566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, gn->parent->gradient));
245cd1c4666STristan Konolige
246cd1c4666STristan Konolige /* Update damping parameter for LM */
247c0b7dd19SHansol Suh if (gn->reg_type == TAOBRGN_REGULARIZATION_LM) {
248cd1c4666STristan Konolige if (iter > 0) {
249cd1c4666STristan Konolige if (gn->fc_old > tao->fc) {
250cd1c4666STristan Konolige gn->lambda = gn->lambda * gn->downhill_lambda_change;
251cd1c4666STristan Konolige } else {
252cd1c4666STristan Konolige /* uphill step */
253cd1c4666STristan Konolige gn->lambda = gn->lambda * gn->uphill_lambda_change;
254cd1c4666STristan Konolige }
255cd1c4666STristan Konolige }
256cd1c4666STristan Konolige gn->fc_old = tao->fc;
257cd1c4666STristan Konolige }
258cd1c4666STristan Konolige
259e1e80dc8SAlp Dener /* Call general purpose update function */
2601baa6e33SBarry Smith if (gn->parent->ops->update) PetscCall((*gn->parent->ops->update)(gn->parent, gn->parent->niter, gn->parent->user_update));
2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
262737f463aSAlp Dener }
263737f463aSAlp Dener
TaoBRGNGetRegularizationType_BRGN(Tao tao,TaoBRGNRegularizationType * type)264c0b7dd19SHansol Suh static PetscErrorCode TaoBRGNGetRegularizationType_BRGN(Tao tao, TaoBRGNRegularizationType *type)
265c0b7dd19SHansol Suh {
266c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
267c0b7dd19SHansol Suh
268c0b7dd19SHansol Suh PetscFunctionBegin;
269c0b7dd19SHansol Suh *type = gn->reg_type;
270c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
271c0b7dd19SHansol Suh }
272c0b7dd19SHansol Suh
273c0b7dd19SHansol Suh /*@
274c0b7dd19SHansol Suh TaoBRGNGetRegularizationType - Get the `TaoBRGNRegularizationType` of a `TAOBRGN`
275c0b7dd19SHansol Suh
276c0b7dd19SHansol Suh Not collective
277c0b7dd19SHansol Suh
278c0b7dd19SHansol Suh Input Parameter:
279c0b7dd19SHansol Suh . tao - a `Tao` of type `TAOBRGN`
280c0b7dd19SHansol Suh
281c0b7dd19SHansol Suh Output Parameter:
282c0b7dd19SHansol Suh . type - the `TaoBRGNRegularizationType`
283c0b7dd19SHansol Suh
284c0b7dd19SHansol Suh Level: advanced
285c0b7dd19SHansol Suh
286c0b7dd19SHansol Suh .seealso: [](ch_tao), `Tao`, `TAOBRGN`, `TaoBRGNRegularizationType`, `TaoBRGNSetRegularizationType()`
287c0b7dd19SHansol Suh @*/
TaoBRGNGetRegularizationType(Tao tao,TaoBRGNRegularizationType * type)288c0b7dd19SHansol Suh PetscErrorCode TaoBRGNGetRegularizationType(Tao tao, TaoBRGNRegularizationType *type)
289c0b7dd19SHansol Suh {
290c0b7dd19SHansol Suh PetscFunctionBegin;
291c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
292c0b7dd19SHansol Suh PetscAssertPointer(type, 2);
293c0b7dd19SHansol Suh PetscUseMethod((PetscObject)tao, "TaoBRGNGetRegularizationType_C", (Tao, TaoBRGNRegularizationType *), (tao, type));
294c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
295c0b7dd19SHansol Suh }
296c0b7dd19SHansol Suh
TaoBRGNSetRegularizationType_BRGN(Tao tao,TaoBRGNRegularizationType type)297c0b7dd19SHansol Suh static PetscErrorCode TaoBRGNSetRegularizationType_BRGN(Tao tao, TaoBRGNRegularizationType type)
298c0b7dd19SHansol Suh {
299c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
300c0b7dd19SHansol Suh
301c0b7dd19SHansol Suh PetscFunctionBegin;
302c0b7dd19SHansol Suh gn->reg_type = type;
303c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
304c0b7dd19SHansol Suh }
305c0b7dd19SHansol Suh
306c0b7dd19SHansol Suh /*@
307c0b7dd19SHansol Suh TaoBRGNSetRegularizationType - Set the `TaoBRGNRegularizationType` of a `TAOBRGN`
308c0b7dd19SHansol Suh
309c0b7dd19SHansol Suh Logically collective
310c0b7dd19SHansol Suh
311c0b7dd19SHansol Suh Input Parameters:
312c0b7dd19SHansol Suh + tao - a `Tao` of type `TAOBRGN`
313c0b7dd19SHansol Suh - type - the `TaoBRGNRegularizationType`
314c0b7dd19SHansol Suh
315c0b7dd19SHansol Suh Level: advanced
316c0b7dd19SHansol Suh
317c0b7dd19SHansol Suh .seealso: [](ch_tao), `Tao`, `TAOBRGN`, `TaoBRGNRegularizationType`, `TaoBRGNGetRegularizationType`
318c0b7dd19SHansol Suh @*/
TaoBRGNSetRegularizationType(Tao tao,TaoBRGNRegularizationType type)319c0b7dd19SHansol Suh PetscErrorCode TaoBRGNSetRegularizationType(Tao tao, TaoBRGNRegularizationType type)
320c0b7dd19SHansol Suh {
321c0b7dd19SHansol Suh PetscFunctionBegin;
322c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
323c0b7dd19SHansol Suh PetscValidLogicalCollectiveEnum(tao, type, 2);
324c0b7dd19SHansol Suh PetscTryMethod((PetscObject)tao, "TaoBRGNSetRegularizationType_C", (Tao, TaoBRGNRegularizationType), (tao, type));
325c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
326c0b7dd19SHansol Suh }
327c0b7dd19SHansol Suh
TaoSolve_BRGN(Tao tao)328d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BRGN(Tao tao)
329d71ae5a4SJacob Faibussowitsch {
330737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data;
331737f463aSAlp Dener
332737f463aSAlp Dener PetscFunctionBegin;
3339566063dSJacob Faibussowitsch PetscCall(TaoSolve(gn->subsolver));
334e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */
335e1e80dc8SAlp Dener tao->nfuncs = gn->subsolver->nfuncs;
336e1e80dc8SAlp Dener tao->ngrads = gn->subsolver->ngrads;
337e1e80dc8SAlp Dener tao->nfuncgrads = gn->subsolver->nfuncgrads;
338e1e80dc8SAlp Dener tao->nhess = gn->subsolver->nhess;
339e1e80dc8SAlp Dener tao->niter = gn->subsolver->niter;
340e1e80dc8SAlp Dener tao->ksp_its = gn->subsolver->ksp_its;
341e1e80dc8SAlp Dener tao->ksp_tot_its = gn->subsolver->ksp_tot_its;
3429566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(gn->subsolver, &tao->reason));
343e1e80dc8SAlp Dener /* Update vectors */
3449566063dSJacob Faibussowitsch PetscCall(VecCopy(gn->subsolver->solution, tao->solution));
3459566063dSJacob Faibussowitsch PetscCall(VecCopy(gn->subsolver->gradient, tao->gradient));
3463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
347737f463aSAlp Dener }
348737f463aSAlp Dener
TaoSetFromOptions_BRGN(Tao tao,PetscOptionItems PetscOptionsObject)349ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BRGN(Tao tao, PetscOptionItems PetscOptionsObject)
350d71ae5a4SJacob Faibussowitsch {
351737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data;
352cd1c4666STristan Konolige TaoLineSearch ls;
353737f463aSAlp Dener
354737f463aSAlp Dener PetscFunctionBegin;
355d0609cedSBarry 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.");
3569566063dSJacob 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));
3579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_brgn_regularizer_weight", "regularizer weight (default 1e-4)", "", gn->lambda, &gn->lambda, NULL));
3589566063dSJacob 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));
3599566063dSJacob 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));
3609566063dSJacob 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));
361c0b7dd19SHansol Suh PetscCall(PetscOptionsEnum("-tao_brgn_regularization_type", "regularization type", "", TaoBRGNRegularizationTypes, (PetscEnum)gn->reg_type, (PetscEnum *)&gn->reg_type, NULL));
362d0609cedSBarry Smith PetscOptionsHeadEnd();
363cd1c4666STristan Konolige /* set unit line search direction as the default when using the lm regularizer */
364c0b7dd19SHansol Suh if (gn->reg_type == TAOBRGN_REGULARIZATION_LM) {
3659566063dSJacob Faibussowitsch PetscCall(TaoGetLineSearch(gn->subsolver, &ls));
3669566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls, TAOLINESEARCHUNIT));
367cd1c4666STristan Konolige }
3689566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(gn->subsolver));
3693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
370737f463aSAlp Dener }
371737f463aSAlp Dener
TaoView_BRGN(Tao tao,PetscViewer viewer)372d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BRGN(Tao tao, PetscViewer viewer)
373d71ae5a4SJacob Faibussowitsch {
374737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data;
375c0b7dd19SHansol Suh PetscBool isascii;
376737f463aSAlp Dener
377737f463aSAlp Dener PetscFunctionBegin;
378c0b7dd19SHansol Suh PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
379c0b7dd19SHansol Suh if (isascii) {
380c0b7dd19SHansol Suh PetscCall(PetscViewerASCIIPushTab(viewer));
381c0b7dd19SHansol Suh PetscCall(PetscViewerASCIIPrintf(viewer, "Regularizer weight: %g\n", (double)gn->lambda));
382c0b7dd19SHansol Suh PetscCall(PetscViewerASCIIPrintf(viewer, "BRGN Regularization Type: %s\n", TaoBRGNRegularizationTypes[gn->reg_type]));
383c0b7dd19SHansol Suh switch (gn->reg_type) {
384c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L1DICT:
385c0b7dd19SHansol Suh PetscCall(PetscViewerASCIIPrintf(viewer, "L1 smooth epsilon: %g\n", (double)gn->epsilon));
386c0b7dd19SHansol Suh break;
387c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_LM:
388c0b7dd19SHansol Suh PetscCall(PetscViewerASCIIPrintf(viewer, "Downhill trust region decrease factor:: %g\n", (double)gn->downhill_lambda_change));
389c0b7dd19SHansol Suh PetscCall(PetscViewerASCIIPrintf(viewer, "Uphill trust region increase factor:: %g\n", (double)gn->uphill_lambda_change));
390c0b7dd19SHansol Suh break;
391c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PROX:
392c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_L2PURE:
393c0b7dd19SHansol Suh case TAOBRGN_REGULARIZATION_USER:
394c0b7dd19SHansol Suh default:
395c0b7dd19SHansol Suh break;
396c0b7dd19SHansol Suh }
397c0b7dd19SHansol Suh PetscCall(PetscViewerASCIIPopTab(viewer));
398c0b7dd19SHansol Suh }
3999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
4009566063dSJacob Faibussowitsch PetscCall(TaoView(gn->subsolver, viewer));
4019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
403737f463aSAlp Dener }
404737f463aSAlp Dener
TaoSetUp_BRGN(Tao tao)405d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BRGN(Tao tao)
406d71ae5a4SJacob Faibussowitsch {
407737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data;
408737f463aSAlp Dener PetscBool is_bnls, is_bntr, is_bntl;
409*a336c150SZach Atkins PetscInt n, N, K; /* dict has size K*N*/
410737f463aSAlp Dener
411737f463aSAlp Dener PetscFunctionBegin;
4123c859ba3SBarry Smith PetscCheck(tao->ls_res, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetResidualRoutine() must be called before setup!");
4139566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNLS, &is_bnls));
4149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNTR, &is_bntr));
4159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)gn->subsolver, TAOBNTL, &is_bntl));
4163c859ba3SBarry Smith PetscCheck((!is_bnls && !is_bntr && !is_bntl) || tao->ls_jac, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetResidualJacobianRoutine() must be called before setup!");
41748a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
41848a46eb9SPierre Jolivet if (!gn->x_work) PetscCall(VecDuplicate(tao->solution, &gn->x_work));
41948a46eb9SPierre Jolivet if (!gn->r_work) PetscCall(VecDuplicate(tao->ls_res, &gn->r_work));
420e1e80dc8SAlp Dener if (!gn->x_old) {
4219566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &gn->x_old));
4229566063dSJacob Faibussowitsch PetscCall(VecSet(gn->x_old, 0.0));
423e1e80dc8SAlp Dener }
4247cea06e1SXiang Huang
425c0b7dd19SHansol Suh if (TAOBRGN_REGULARIZATION_L1DICT == gn->reg_type) {
4262036730cSSajid Ali if (!gn->y) {
42730eeff36SXiang Huang if (gn->D) {
4289566063dSJacob 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 */
4299566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(gn->D, NULL, &gn->y));
43030eeff36SXiang Huang } else {
431aaa8cc7dSPierre Jolivet PetscCall(VecDuplicate(tao->solution, &gn->y)); /* If user does not setup dict matrix, use identity matrix, K=N */
43230eeff36SXiang Huang }
4339566063dSJacob Faibussowitsch PetscCall(VecSet(gn->y, 0.0));
4347cea06e1SXiang Huang }
43548a46eb9SPierre Jolivet if (!gn->y_work) PetscCall(VecDuplicate(gn->y, &gn->y_work));
4368ac80d48SXiang Huang if (!gn->diag) {
4379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(gn->y, &gn->diag));
4389566063dSJacob Faibussowitsch PetscCall(VecSet(gn->diag, 0.0));
4398ac80d48SXiang Huang }
44030eeff36SXiang Huang }
441c0b7dd19SHansol Suh if (TAOBRGN_REGULARIZATION_LM == gn->reg_type) {
44248a46eb9SPierre Jolivet if (!gn->diag) PetscCall(MatCreateVecs(tao->ls_jac, &gn->diag, NULL));
44348a46eb9SPierre Jolivet if (!gn->damping) PetscCall(MatCreateVecs(tao->ls_jac, &gn->damping, NULL));
444cd1c4666STristan Konolige }
4450d71dc2bSXiang Huang
446e1e80dc8SAlp Dener if (!tao->setupcalled) {
447737f463aSAlp Dener /* Hessian setup */
4485eb5f4d6SAlp Dener if (gn->mat_explicit) {
4499566063dSJacob Faibussowitsch PetscCall(TaoComputeResidualJacobian(tao, tao->solution, tao->ls_jac, tao->ls_jac_pre));
450fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(tao->ls_jac, tao->ls_jac, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &gn->H));
4515eb5f4d6SAlp Dener } else {
4529566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n));
4539566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N));
4549566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &gn->H));
4559566063dSJacob Faibussowitsch PetscCall(MatSetSizes(gn->H, n, n, N, N));
4569566063dSJacob Faibussowitsch PetscCall(MatSetType(gn->H, MATSHELL));
4579566063dSJacob Faibussowitsch PetscCall(MatSetOption(gn->H, MAT_SYMMETRIC, PETSC_TRUE));
45857d50842SBarry Smith PetscCall(MatShellSetOperation(gn->H, MATOP_MULT, (PetscErrorCodeFn *)GNHessianProd));
4599566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(gn->H, gn));
4605eb5f4d6SAlp Dener }
4619566063dSJacob Faibussowitsch PetscCall(MatSetUp(gn->H));
462a5b23f4aSJose E. Roman /* Subsolver setup,include initial vector and dictionary D */
4639566063dSJacob Faibussowitsch PetscCall(TaoSetUpdate(gn->subsolver, GNHookFunction, gn));
4649566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(gn->subsolver, tao->solution));
4651baa6e33SBarry Smith if (tao->bounded) PetscCall(TaoSetVariableBounds(gn->subsolver, tao->XL, tao->XU));
4669566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(gn->subsolver, tao->ls_res, tao->ops->computeresidual, tao->user_lsresP));
4679566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianResidualRoutine(gn->subsolver, tao->ls_jac, tao->ls_jac, tao->ops->computeresidualjacobian, tao->user_lsjacP));
4689566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(gn->subsolver, NULL, GNObjectiveGradientEval, gn));
4699566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(gn->subsolver, gn->H, gn->H, GNComputeHessian, gn));
470e1e80dc8SAlp Dener /* Propagate some options down */
4719566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(gn->subsolver, tao->gatol, tao->grtol, tao->gttol));
4729566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(gn->subsolver, tao->max_it));
4739566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumFunctionEvaluations(gn->subsolver, tao->max_funcs));
4749566063dSJacob Faibussowitsch PetscCall(TaoSetUp(gn->subsolver));
475e1e80dc8SAlp Dener }
4763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
477737f463aSAlp Dener }
478737f463aSAlp Dener
TaoDestroy_BRGN(Tao tao)479d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BRGN(Tao tao)
480d71ae5a4SJacob Faibussowitsch {
481737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data;
482737f463aSAlp Dener
483737f463aSAlp Dener PetscFunctionBegin;
484737f463aSAlp Dener if (tao->setupcalled) {
4859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->gradient));
4869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->x_work));
4879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->r_work));
4889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->x_old));
4899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->diag));
4909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->y));
4919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->y_work));
492737f463aSAlp Dener }
4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->damping));
4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gn->diag));
4959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&gn->H));
4969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&gn->D));
4979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&gn->Hreg));
4989566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&gn->subsolver));
499e1e80dc8SAlp Dener gn->parent = NULL;
5009566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
501c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNGetRegularizationType_C", NULL));
502c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizationType_C", NULL));
503c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNGetDampingVector_C", NULL));
504c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetDictionaryMatrix_C", NULL));
505c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNGetSubsolver_C", NULL));
506c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizerWeight_C", NULL));
507c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetL1SmoothEpsilon_C", NULL));
508c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizerObjectiveAndGradientRoutine_C", NULL));
509c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizerHessianRoutine_C", NULL));
510c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
511c0b7dd19SHansol Suh }
512c0b7dd19SHansol Suh
513c0b7dd19SHansol Suh /*@
514c0b7dd19SHansol Suh TaoBRGNGetSubsolver - Get the pointer to the subsolver inside a `TAOBRGN`
515c0b7dd19SHansol Suh
516c0b7dd19SHansol Suh Collective
517c0b7dd19SHansol Suh
518c0b7dd19SHansol Suh Input Parameters:
519c0b7dd19SHansol Suh + tao - the Tao solver context
520c0b7dd19SHansol Suh - subsolver - the `Tao` sub-solver context
521c0b7dd19SHansol Suh
522c0b7dd19SHansol Suh Level: advanced
523c0b7dd19SHansol Suh
524c0b7dd19SHansol Suh .seealso: `Tao`, `Mat`, `TAOBRGN`
525c0b7dd19SHansol Suh @*/
TaoBRGNGetSubsolver(Tao tao,Tao * subsolver)526c0b7dd19SHansol Suh PetscErrorCode TaoBRGNGetSubsolver(Tao tao, Tao *subsolver)
527c0b7dd19SHansol Suh {
528c0b7dd19SHansol Suh PetscFunctionBegin;
529c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
530c0b7dd19SHansol Suh PetscUseMethod((PetscObject)tao, "TaoBRGNGetSubsolver_C", (Tao, Tao *), (tao, subsolver));
531c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
532c0b7dd19SHansol Suh }
533c0b7dd19SHansol Suh
TaoBRGNGetSubsolver_BRGN(Tao tao,Tao * subsolver)534c0b7dd19SHansol Suh static PetscErrorCode TaoBRGNGetSubsolver_BRGN(Tao tao, Tao *subsolver)
535c0b7dd19SHansol Suh {
536c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
537c0b7dd19SHansol Suh
538c0b7dd19SHansol Suh PetscFunctionBegin;
539c0b7dd19SHansol Suh *subsolver = gn->subsolver;
540c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
541c0b7dd19SHansol Suh }
542c0b7dd19SHansol Suh
543c0b7dd19SHansol Suh /*@
544c0b7dd19SHansol Suh TaoBRGNSetRegularizerWeight - Set the regularizer weight for the Gauss-Newton least-squares algorithm
545c0b7dd19SHansol Suh
546c0b7dd19SHansol Suh Collective
547c0b7dd19SHansol Suh
548c0b7dd19SHansol Suh Input Parameters:
549c0b7dd19SHansol Suh + tao - the `Tao` solver context
550c0b7dd19SHansol Suh - lambda - L1-norm regularizer weight
551c0b7dd19SHansol Suh
552c0b7dd19SHansol Suh Level: beginner
553c0b7dd19SHansol Suh
554c0b7dd19SHansol Suh .seealso: `Tao`, `Mat`, `TAOBRGN`
555c0b7dd19SHansol Suh @*/
TaoBRGNSetRegularizerWeight(Tao tao,PetscReal lambda)556c0b7dd19SHansol Suh PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao, PetscReal lambda)
557c0b7dd19SHansol Suh {
558c0b7dd19SHansol Suh PetscFunctionBegin;
559c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
560c0b7dd19SHansol Suh PetscValidLogicalCollectiveReal(tao, lambda, 2);
561c0b7dd19SHansol Suh PetscTryMethod((PetscObject)tao, "TaoBRGNSetRegularizerWeight_C", (Tao, PetscReal), (tao, lambda));
562c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
563c0b7dd19SHansol Suh }
564c0b7dd19SHansol Suh
TaoBRGNSetRegularizerWeight_BRGN(Tao tao,PetscReal lambda)565c0b7dd19SHansol Suh static PetscErrorCode TaoBRGNSetRegularizerWeight_BRGN(Tao tao, PetscReal lambda)
566c0b7dd19SHansol Suh {
567c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
568c0b7dd19SHansol Suh
569c0b7dd19SHansol Suh PetscFunctionBegin;
570c0b7dd19SHansol Suh gn->lambda = lambda;
571c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
572c0b7dd19SHansol Suh }
573c0b7dd19SHansol Suh
574c0b7dd19SHansol Suh /*@
575c0b7dd19SHansol Suh TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm
576c0b7dd19SHansol Suh
577c0b7dd19SHansol Suh Collective
578c0b7dd19SHansol Suh
579c0b7dd19SHansol Suh Input Parameters:
580c0b7dd19SHansol Suh + tao - the `Tao` solver context
581c0b7dd19SHansol Suh - epsilon - L1-norm smooth approximation parameter
582c0b7dd19SHansol Suh
583c0b7dd19SHansol Suh Level: advanced
584c0b7dd19SHansol Suh
585c0b7dd19SHansol Suh .seealso: `Tao`, `Mat`, `TAOBRGN`
586c0b7dd19SHansol Suh @*/
TaoBRGNSetL1SmoothEpsilon(Tao tao,PetscReal epsilon)587c0b7dd19SHansol Suh PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao, PetscReal epsilon)
588c0b7dd19SHansol Suh {
589c0b7dd19SHansol Suh PetscFunctionBegin;
590c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
591c0b7dd19SHansol Suh PetscValidLogicalCollectiveReal(tao, epsilon, 2);
592c0b7dd19SHansol Suh PetscTryMethod((PetscObject)tao, "TaoBRGNSetL1SmoothEpsilon_C", (Tao, PetscReal), (tao, epsilon));
593c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
594c0b7dd19SHansol Suh }
595c0b7dd19SHansol Suh
TaoBRGNSetL1SmoothEpsilon_BRGN(Tao tao,PetscReal epsilon)596c0b7dd19SHansol Suh static PetscErrorCode TaoBRGNSetL1SmoothEpsilon_BRGN(Tao tao, PetscReal epsilon)
597c0b7dd19SHansol Suh {
598c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
599c0b7dd19SHansol Suh
600c0b7dd19SHansol Suh PetscFunctionBegin;
601c0b7dd19SHansol Suh gn->epsilon = epsilon;
602c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
603c0b7dd19SHansol Suh }
604c0b7dd19SHansol Suh
605c0b7dd19SHansol Suh /*@
606c0b7dd19SHansol Suh TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem)
607c0b7dd19SHansol Suh
608c0b7dd19SHansol Suh Input Parameters:
609c0b7dd19SHansol Suh + tao - the `Tao` context
610c0b7dd19SHansol Suh - dict - the user specified dictionary matrix. We allow to set a `NULL` dictionary, which means identity matrix by default
611c0b7dd19SHansol Suh
612c0b7dd19SHansol Suh Level: advanced
613c0b7dd19SHansol Suh
614c0b7dd19SHansol Suh .seealso: `Tao`, `Mat`, `TAOBRGN`
615c0b7dd19SHansol Suh @*/
TaoBRGNSetDictionaryMatrix(Tao tao,Mat dict)616c0b7dd19SHansol Suh PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao, Mat dict)
617c0b7dd19SHansol Suh {
618c0b7dd19SHansol Suh PetscFunctionBegin;
619c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
620c0b7dd19SHansol Suh PetscTryMethod((PetscObject)tao, "TaoBRGNSetDictionaryMatrix_C", (Tao, Mat), (tao, dict));
621c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
622c0b7dd19SHansol Suh }
623c0b7dd19SHansol Suh
TaoBRGNSetDictionaryMatrix_BRGN(Tao tao,Mat dict)624c0b7dd19SHansol Suh static PetscErrorCode TaoBRGNSetDictionaryMatrix_BRGN(Tao tao, Mat dict)
625c0b7dd19SHansol Suh {
626c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
627c0b7dd19SHansol Suh
628c0b7dd19SHansol Suh PetscFunctionBegin;
629c0b7dd19SHansol Suh if (dict) {
630c0b7dd19SHansol Suh PetscValidHeaderSpecific(dict, MAT_CLASSID, 2);
631c0b7dd19SHansol Suh PetscCheckSameComm(tao, 1, dict, 2);
632c0b7dd19SHansol Suh PetscCall(PetscObjectReference((PetscObject)dict));
633c0b7dd19SHansol Suh }
634c0b7dd19SHansol Suh PetscCall(MatDestroy(&gn->D));
635c0b7dd19SHansol Suh gn->D = dict;
636c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
637c0b7dd19SHansol Suh }
638c0b7dd19SHansol Suh
639c0b7dd19SHansol Suh /*@C
640c0b7dd19SHansol Suh TaoBRGNSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back
641c0b7dd19SHansol Suh function into the algorithm.
642c0b7dd19SHansol Suh
643c0b7dd19SHansol Suh Input Parameters:
644c0b7dd19SHansol Suh + tao - the Tao context
645c0b7dd19SHansol Suh . func - function pointer for the regularizer value and gradient evaluation
646c0b7dd19SHansol Suh - ctx - user context for the regularizer
647c0b7dd19SHansol Suh
648c0b7dd19SHansol Suh Calling sequence:
649c0b7dd19SHansol Suh + tao - the `Tao` context
650c0b7dd19SHansol Suh . u - the location at which to compute the objective and gradient
651c0b7dd19SHansol Suh . val - location to store objective function value
652c0b7dd19SHansol Suh . g - location to store gradient
653c0b7dd19SHansol Suh - ctx - user context for the regularizer Hessian
654c0b7dd19SHansol Suh
655c0b7dd19SHansol Suh Level: advanced
656c0b7dd19SHansol Suh
657c0b7dd19SHansol Suh .seealso: `Tao`, `Mat`, `TAOBRGN`
658c0b7dd19SHansol Suh @*/
TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (* func)(Tao tao,Vec u,PetscReal * val,Vec g,PetscCtx ctx),PetscCtx ctx)6592a8381b2SBarry Smith PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao tao, Vec u, PetscReal *val, Vec g, PetscCtx ctx), PetscCtx ctx)
660c0b7dd19SHansol Suh {
661c0b7dd19SHansol Suh PetscFunctionBegin;
662c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
663c0b7dd19SHansol Suh PetscTryMethod((PetscObject)tao, "TaoBRGNSetRegularizerObjectiveAndGradientRoutine_C", (Tao, PetscErrorCode (*)(Tao, Vec, PetscReal *, Vec, void *), void *), (tao, func, ctx));
664c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
665c0b7dd19SHansol Suh }
666c0b7dd19SHansol Suh
TaoBRGNSetRegularizerObjectiveAndGradientRoutine_BRGN(Tao tao,PetscErrorCode (* func)(Tao tao,Vec u,PetscReal * val,Vec g,PetscCtx ctx),PetscCtx ctx)6672a8381b2SBarry Smith static PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine_BRGN(Tao tao, PetscErrorCode (*func)(Tao tao, Vec u, PetscReal *val, Vec g, PetscCtx ctx), PetscCtx ctx)
668c0b7dd19SHansol Suh {
669c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
670c0b7dd19SHansol Suh
671c0b7dd19SHansol Suh PetscFunctionBegin;
672c0b7dd19SHansol Suh if (ctx) gn->reg_obj_ctx = ctx;
673c0b7dd19SHansol Suh if (func) gn->regularizerobjandgrad = func;
674c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
675c0b7dd19SHansol Suh }
676c0b7dd19SHansol Suh
677c0b7dd19SHansol Suh /*@C
678c0b7dd19SHansol Suh TaoBRGNSetRegularizerHessianRoutine - Sets the user-defined regularizer call-back
679c0b7dd19SHansol Suh function into the algorithm.
680c0b7dd19SHansol Suh
681c0b7dd19SHansol Suh Input Parameters:
682c0b7dd19SHansol Suh + tao - the `Tao` context
683c0b7dd19SHansol Suh . Hreg - user-created matrix for the Hessian of the regularization term
684c0b7dd19SHansol Suh . func - function pointer for the regularizer Hessian evaluation
685c0b7dd19SHansol Suh - ctx - user context for the regularizer Hessian
686c0b7dd19SHansol Suh
687c0b7dd19SHansol Suh Calling sequence:
688c0b7dd19SHansol Suh + tao - the `Tao` context
689c0b7dd19SHansol Suh . u - the location at which to compute the Hessian
690c0b7dd19SHansol Suh . Hreg - user-created matrix for the Hessian of the regularization term
691c0b7dd19SHansol Suh - ctx - user context for the regularizer Hessian
692c0b7dd19SHansol Suh
693c0b7dd19SHansol Suh Level: advanced
694c0b7dd19SHansol Suh
695c0b7dd19SHansol Suh .seealso: `Tao`, `Mat`, `TAOBRGN`
696c0b7dd19SHansol Suh @*/
TaoBRGNSetRegularizerHessianRoutine(Tao tao,Mat Hreg,PetscErrorCode (* func)(Tao tao,Vec u,Mat Hreg,PetscCtx ctx),PetscCtx ctx)6972a8381b2SBarry Smith PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao, Mat Hreg, PetscErrorCode (*func)(Tao tao, Vec u, Mat Hreg, PetscCtx ctx), PetscCtx ctx)
698c0b7dd19SHansol Suh {
699c0b7dd19SHansol Suh PetscFunctionBegin;
700c0b7dd19SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
701c0b7dd19SHansol Suh PetscTryMethod((PetscObject)tao, "TaoBRGNSetRegularizerHessianRoutine_C", (Tao, Mat, PetscErrorCode (*)(Tao, Vec, Mat, void *), void *), (tao, Hreg, func, ctx));
702c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
703c0b7dd19SHansol Suh }
704c0b7dd19SHansol Suh
TaoBRGNSetRegularizerHessianRoutine_BRGN(Tao tao,Mat Hreg,PetscErrorCode (* func)(Tao tao,Vec u,Mat Hreg,PetscCtx ctx),PetscCtx ctx)7052a8381b2SBarry Smith static PetscErrorCode TaoBRGNSetRegularizerHessianRoutine_BRGN(Tao tao, Mat Hreg, PetscErrorCode (*func)(Tao tao, Vec u, Mat Hreg, PetscCtx ctx), PetscCtx ctx)
706c0b7dd19SHansol Suh {
707c0b7dd19SHansol Suh TAO_BRGN *gn = (TAO_BRGN *)tao->data;
708c0b7dd19SHansol Suh
709c0b7dd19SHansol Suh PetscFunctionBegin;
710c0b7dd19SHansol Suh if (Hreg) {
711c0b7dd19SHansol Suh PetscValidHeaderSpecific(Hreg, MAT_CLASSID, 2);
712c0b7dd19SHansol Suh PetscCheckSameComm(tao, 1, Hreg, 2);
713c0b7dd19SHansol Suh } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "NULL Hessian detected! User must provide valid Hessian for the regularizer.");
714c0b7dd19SHansol Suh if (ctx) gn->reg_hess_ctx = ctx;
715c0b7dd19SHansol Suh if (func) gn->regularizerhessian = func;
716c0b7dd19SHansol Suh if (Hreg) {
717c0b7dd19SHansol Suh PetscCall(PetscObjectReference((PetscObject)Hreg));
718c0b7dd19SHansol Suh PetscCall(MatDestroy(&gn->Hreg));
719c0b7dd19SHansol Suh gn->Hreg = Hreg;
720c0b7dd19SHansol Suh }
7213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
722737f463aSAlp Dener }
723737f463aSAlp Dener
7243850be85SAlp Dener /*MC
7253850be85SAlp Dener TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares
7262fe279fdSBarry Smith problems with bound constraints. This algorithm is a thin wrapper around `TAOBNTL`
727463fc0ecSAlp Dener that constructs the Gauss-Newton problem with the user-provided least-squares
72860bb7533SHansol Suh residual and Jacobian. The algorithm offers an L2-norm ("l2pure"), L2-norm proximal point ("l2prox")
72960bb7533SHansol Suh regularizer, and L1-norm dictionary regularizer ("l1dict"), where we approximate the
73001b716f5SXiang Huang L1-norm ||x||_1 by sum_i(sqrt(x_i^2+epsilon^2)-epsilon) with a small positive number epsilon.
731cd1c4666STristan Konolige Also offered is the "lm" regularizer which uses a scaled diagonal of J^T J.
7322fe279fdSBarry Smith With the "lm" regularizer, `TAOBRGN` is a Levenberg-Marquardt optimizer.
73301b716f5SXiang Huang The user can also provide own regularization function.
7343850be85SAlp Dener
7353850be85SAlp Dener Options Database Keys:
736cd1c4666STristan Konolige + -tao_brgn_regularization_type - regularization type ("user", "l2prox", "l2pure", "l1dict", "lm") (default "l2prox")
737c061e8e2SXiang Huang . -tao_brgn_regularizer_weight - regularizer weight (default 1e-4)
738c061e8e2SXiang Huang - -tao_brgn_l1_smooth_epsilon - L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6)
7393850be85SAlp Dener
7403850be85SAlp Dener Level: beginner
7412fe279fdSBarry Smith
7422fe279fdSBarry Smith .seealso: `Tao`, `TaoBRGNGetSubsolver()`, `TaoBRGNSetRegularizerWeight()`, `TaoBRGNSetL1SmoothEpsilon()`, `TaoBRGNSetDictionaryMatrix()`,
7432fe279fdSBarry Smith `TaoBRGNSetRegularizerObjectiveAndGradientRoutine()`, `TaoBRGNSetRegularizerHessianRoutine()`
7443850be85SAlp Dener M*/
TaoCreate_BRGN(Tao tao)745d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao)
746d71ae5a4SJacob Faibussowitsch {
747737f463aSAlp Dener TAO_BRGN *gn;
748737f463aSAlp Dener
749737f463aSAlp Dener PetscFunctionBegin;
7504dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&gn));
751737f463aSAlp Dener
752737f463aSAlp Dener tao->ops->destroy = TaoDestroy_BRGN;
753737f463aSAlp Dener tao->ops->setup = TaoSetUp_BRGN;
754737f463aSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BRGN;
755737f463aSAlp Dener tao->ops->view = TaoView_BRGN;
756737f463aSAlp Dener tao->ops->solve = TaoSolve_BRGN;
757737f463aSAlp Dener
758606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
759606f75f6SBarry Smith
7603ec1f749SStefano Zampini tao->data = gn;
761c0b7dd19SHansol Suh gn->reg_type = TAOBRGN_REGULARIZATION_L2PROX;
762e1e80dc8SAlp Dener gn->lambda = 1e-4;
7638ac80d48SXiang Huang gn->epsilon = 1e-6;
764cd1c4666STristan Konolige gn->downhill_lambda_change = 1. / 5.;
765cd1c4666STristan Konolige gn->uphill_lambda_change = 1.5;
766e1e80dc8SAlp Dener gn->parent = tao;
767737f463aSAlp Dener
7689566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &gn->subsolver));
7699566063dSJacob Faibussowitsch PetscCall(TaoSetType(gn->subsolver, TAOBNLS));
7709566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(gn->subsolver, "tao_brgn_subsolver_"));
771c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNGetRegularizationType_C", TaoBRGNGetRegularizationType_BRGN));
772c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizationType_C", TaoBRGNSetRegularizationType_BRGN));
773c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNGetDampingVector_C", TaoBRGNGetDampingVector_BRGN));
774c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetDictionaryMatrix_C", TaoBRGNSetDictionaryMatrix_BRGN));
775c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNGetSubsolver_C", TaoBRGNGetSubsolver_BRGN));
776c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizerWeight_C", TaoBRGNSetRegularizerWeight_BRGN));
777c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetL1SmoothEpsilon_C", TaoBRGNSetL1SmoothEpsilon_BRGN));
778c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizerObjectiveAndGradientRoutine_C", TaoBRGNSetRegularizerObjectiveAndGradientRoutine_BRGN));
779c0b7dd19SHansol Suh PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoBRGNSetRegularizerHessianRoutine_C", TaoBRGNSetRegularizerHessianRoutine_BRGN));
7803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
781a3c390cfSAlp Dener }
782