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