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