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