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