1463fc0ecSAlp Dener #include <../src/tao/leastsquares/impls/brgn/brgn.h> /*I "petsctao.h" I*/ 2737f463aSAlp Dener 3470ec3f8SXiang Huang #define BRGN_REGULARIZATION_USER 0 4470ec3f8SXiang Huang #define BRGN_REGULARIZATION_L2PROX 1 5a1c74439SHansol Suh #define BRGN_REGULARIZATION_L2PURE 2 6a1c74439SHansol Suh #define BRGN_REGULARIZATION_L1DICT 3 7*cd1c4666STristan Konolige #define BRGN_REGULARIZATION_LM 4 8*cd1c4666STristan Konolige #define BRGN_REGULARIZATION_TYPES 5 9a3c390cfSAlp Dener 10*cd1c4666STristan Konolige static const char *BRGN_REGULARIZATION_TABLE[64] = {"user","l2prox","l2pure","l1dict","lm"}; 11a3c390cfSAlp Dener 120d71dc2bSXiang Huang static PetscErrorCode GNHessianProd(Mat H,Vec in,Vec out) 130d71dc2bSXiang Huang { 140d71dc2bSXiang Huang TAO_BRGN *gn; 150d71dc2bSXiang Huang PetscErrorCode ierr; 160d71dc2bSXiang Huang 170d71dc2bSXiang Huang PetscFunctionBegin; 180d71dc2bSXiang Huang ierr = MatShellGetContext(H,&gn);CHKERRQ(ierr); 190d71dc2bSXiang Huang ierr = MatMult(gn->subsolver->ls_jac,in,gn->r_work);CHKERRQ(ierr); 200d71dc2bSXiang Huang ierr = MatMultTranspose(gn->subsolver->ls_jac,gn->r_work,out);CHKERRQ(ierr); 21a3c390cfSAlp Dener switch (gn->reg_type) { 22470ec3f8SXiang Huang case BRGN_REGULARIZATION_USER: 23a3c390cfSAlp Dener ierr = MatMult(gn->Hreg,in,gn->x_work);CHKERRQ(ierr); 248ac80d48SXiang Huang ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr); 25a3c390cfSAlp Dener break; 26a1c74439SHansol Suh case BRGN_REGULARIZATION_L2PURE: 27a1c74439SHansol Suh ierr = VecAXPY(out,gn->lambda,in);CHKERRQ(ierr); 28a1c74439SHansol Suh break; 29470ec3f8SXiang Huang case BRGN_REGULARIZATION_L2PROX: 30a3c390cfSAlp Dener ierr = VecAXPY(out,gn->lambda,in);CHKERRQ(ierr); 31a3c390cfSAlp Dener break; 32470ec3f8SXiang Huang case BRGN_REGULARIZATION_L1DICT: 33a3c390cfSAlp Dener /* out = out + lambda*D'*(diag.*(D*in)) */ 34a3c390cfSAlp Dener if (gn->D) { 35a3c390cfSAlp Dener ierr = MatMult(gn->D,in,gn->y);CHKERRQ(ierr);/* y = D*in */ 36a3c390cfSAlp Dener } else { 37a3c390cfSAlp Dener ierr = VecCopy(in,gn->y);CHKERRQ(ierr); 38a3c390cfSAlp Dener } 39a3c390cfSAlp Dener ierr = VecPointwiseMult(gn->y_work,gn->diag,gn->y);CHKERRQ(ierr); /* y_work = diag.*(D*in), where diag = epsilon^2 ./ sqrt(x.^2+epsilon^2).^3 */ 40a3c390cfSAlp Dener if (gn->D) { 41a3c390cfSAlp Dener ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr); /* x_work = D'*(diag.*(D*in)) */ 42a3c390cfSAlp Dener } else { 43a3c390cfSAlp Dener ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr); 44a3c390cfSAlp Dener } 45a3c390cfSAlp Dener ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr); 46a3c390cfSAlp Dener break; 47*cd1c4666STristan Konolige case BRGN_REGULARIZATION_LM: 48*cd1c4666STristan Konolige ierr = VecPointwiseMult(gn->x_work,gn->damping,in);CHKERRQ(ierr); 49*cd1c4666STristan Konolige ierr = VecAXPY(out,1,gn->x_work);CHKERRQ(ierr); 50*cd1c4666STristan Konolige break; 51a3c390cfSAlp Dener } 520d71dc2bSXiang Huang PetscFunctionReturn(0); 530d71dc2bSXiang Huang } 54*cd1c4666STristan Konolige static PetscErrorCode ComputeDamping(TAO_BRGN *gn) 55*cd1c4666STristan Konolige { 56*cd1c4666STristan Konolige const PetscScalar *diag_ary; 57*cd1c4666STristan Konolige PetscScalar *damping_ary; 58*cd1c4666STristan Konolige PetscInt i,n; 59*cd1c4666STristan Konolige PetscErrorCode ierr; 60*cd1c4666STristan Konolige 61*cd1c4666STristan Konolige PetscFunctionBegin; 62*cd1c4666STristan Konolige /* update damping */ 63*cd1c4666STristan Konolige ierr = VecGetArray(gn->damping,&damping_ary);CHKERRQ(ierr); 64*cd1c4666STristan Konolige ierr = VecGetArrayRead(gn->diag,&diag_ary);CHKERRQ(ierr); 65*cd1c4666STristan Konolige ierr = VecGetLocalSize(gn->damping,&n);CHKERRQ(ierr); 66*cd1c4666STristan Konolige for (i=0; i<n; i++) { 67*cd1c4666STristan Konolige damping_ary[i] = PetscClipInterval(diag_ary[i],PETSC_SQRT_MACHINE_EPSILON,PetscSqrtReal(PETSC_MAX_REAL)); 68*cd1c4666STristan Konolige } 69*cd1c4666STristan Konolige ierr = VecScale(gn->damping,gn->lambda);CHKERRQ(ierr); 70*cd1c4666STristan Konolige ierr = VecRestoreArray(gn->damping,&damping_ary);CHKERRQ(ierr); 71*cd1c4666STristan Konolige ierr = VecRestoreArrayRead(gn->diag,&diag_ary);CHKERRQ(ierr); 72*cd1c4666STristan Konolige PetscFunctionReturn(0); 73*cd1c4666STristan Konolige } 74*cd1c4666STristan Konolige 75*cd1c4666STristan Konolige PetscErrorCode TaoBRGNGetDampingVector(Tao tao,Vec *d) 76*cd1c4666STristan Konolige { 77*cd1c4666STristan Konolige TAO_BRGN *gn = (TAO_BRGN *)tao->data; 78*cd1c4666STristan Konolige 79*cd1c4666STristan Konolige PetscFunctionBegin; 80*cd1c4666STristan Konolige if (gn->reg_type != BRGN_REGULARIZATION_LM) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Damping vector is only available if regularization type is lm."); 81*cd1c4666STristan Konolige *d = gn->damping; 82*cd1c4666STristan Konolige PetscFunctionReturn(0); 83*cd1c4666STristan Konolige } 840d71dc2bSXiang Huang 850d71dc2bSXiang Huang static PetscErrorCode GNObjectiveGradientEval(Tao tao,Vec X,PetscReal *fcn,Vec G,void *ptr) 860d71dc2bSXiang Huang { 870d71dc2bSXiang Huang TAO_BRGN *gn = (TAO_BRGN *)ptr; 888e85b1b3SXiang Huang PetscInt K; /* dimension of D*X */ 897cea06e1SXiang Huang PetscScalar yESum; 900d71dc2bSXiang Huang PetscErrorCode ierr; 91a3c390cfSAlp Dener PetscReal f_reg; 920d71dc2bSXiang Huang 930d71dc2bSXiang Huang PetscFunctionBegin; 948e85b1b3SXiang Huang /* compute objective *fcn*/ 95a3c390cfSAlp Dener /* compute first term 0.5*||ls_res||_2^2 */ 960d71dc2bSXiang Huang ierr = TaoComputeResidual(tao,X,tao->ls_res);CHKERRQ(ierr); 97a3c390cfSAlp Dener ierr = VecDot(tao->ls_res,tao->ls_res,fcn);CHKERRQ(ierr); 98a3c390cfSAlp Dener *fcn *= 0.5; 99a3c390cfSAlp Dener /* compute gradient G */ 100a3c390cfSAlp Dener ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr); 101a3c390cfSAlp Dener ierr = MatMultTranspose(tao->ls_jac,tao->ls_res,G);CHKERRQ(ierr); 102a3c390cfSAlp Dener /* add the regularization contribution */ 103a3c390cfSAlp Dener switch (gn->reg_type) { 104470ec3f8SXiang Huang case BRGN_REGULARIZATION_USER: 105a3c390cfSAlp Dener ierr = (*gn->regularizerobjandgrad)(tao,X,&f_reg,gn->x_work,gn->reg_obj_ctx);CHKERRQ(ierr); 106a3c390cfSAlp Dener *fcn += gn->lambda*f_reg; 107a3c390cfSAlp Dener ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr); 108a3c390cfSAlp Dener break; 109a1c74439SHansol Suh case BRGN_REGULARIZATION_L2PURE: 110a1c74439SHansol Suh /* compute f = f + lambda*0.5*xk'*xk */ 111a1c74439SHansol Suh ierr = VecDot(X,X,&f_reg);CHKERRQ(ierr); 112a1c74439SHansol Suh *fcn += gn->lambda*0.5*f_reg; 113a1c74439SHansol Suh /* compute G = G + lambda*xk */ 114a1c74439SHansol Suh ierr = VecAXPY(G,gn->lambda,X);CHKERRQ(ierr); 115a1c74439SHansol Suh break; 116470ec3f8SXiang Huang case BRGN_REGULARIZATION_L2PROX: 1171fc140a9SXiang Huang /* compute f = f + lambda*0.5*(xk - xkm1)'*(xk - xkm1) */ 1181fc140a9SXiang Huang ierr = VecAXPBYPCZ(gn->x_work,1.0,-1.0,0.0,X,gn->x_old);CHKERRQ(ierr); 119a3c390cfSAlp Dener ierr = VecDot(gn->x_work,gn->x_work,&f_reg);CHKERRQ(ierr); 120a3c390cfSAlp Dener *fcn += gn->lambda*0.5*f_reg; 121a3c390cfSAlp Dener /* compute G = G + lambda*(xk - xkm1) */ 122a3c390cfSAlp Dener ierr = VecAXPBYPCZ(G,gn->lambda,-gn->lambda,1.0,X,gn->x_old);CHKERRQ(ierr); 123a3c390cfSAlp Dener break; 124470ec3f8SXiang Huang case BRGN_REGULARIZATION_L1DICT: 125a3c390cfSAlp Dener /* compute f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x*/ 126a3c390cfSAlp Dener if (gn->D) { 1277cea06e1SXiang Huang ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */ 128a3c390cfSAlp Dener } else { 129a3c390cfSAlp Dener ierr = VecCopy(X,gn->y);CHKERRQ(ierr); 130a3c390cfSAlp Dener } 1317cea06e1SXiang Huang ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr); 1327cea06e1SXiang Huang ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 1337cea06e1SXiang Huang ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr); /* gn->y_work = sqrt(y.^2+epsilon^2) */ 134c3b366b1Sprj- ierr = VecSum(gn->y_work,&yESum);CHKERRQ(ierr); 1358e85b1b3SXiang Huang ierr = VecGetSize(gn->y,&K);CHKERRQ(ierr); 136a3c390cfSAlp Dener *fcn += gn->lambda*(yESum - K*gn->epsilon); 1377cea06e1SXiang Huang /* compute G = G + lambda*D'*(y./sqrt(y.^2+epsilon^2)),where y = D*x */ 1387cea06e1SXiang Huang ierr = VecPointwiseDivide(gn->y_work,gn->y,gn->y_work);CHKERRQ(ierr); /* reuse y_work = y./sqrt(y.^2+epsilon^2) */ 139a3c390cfSAlp Dener if (gn->D) { 1407cea06e1SXiang Huang ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr); 141a3c390cfSAlp Dener } else { 142a3c390cfSAlp Dener ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr); 143a3c390cfSAlp Dener } 1448ac80d48SXiang Huang ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr); 145a3c390cfSAlp Dener break; 146a3c390cfSAlp Dener } 1470d71dc2bSXiang Huang PetscFunctionReturn(0); 1480d71dc2bSXiang Huang } 1490d71dc2bSXiang Huang 150737f463aSAlp Dener static PetscErrorCode GNComputeHessian(Tao tao,Vec X,Mat H,Mat Hpre,void *ptr) 151737f463aSAlp Dener { 1528ac80d48SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)ptr; 153*cd1c4666STristan Konolige PetscInt i,n,cstart,cend; 154*cd1c4666STristan Konolige PetscScalar *cnorms,*diag_ary; 155737f463aSAlp Dener PetscErrorCode ierr; 156737f463aSAlp Dener 157737f463aSAlp Dener PetscFunctionBegin; 158e1e80dc8SAlp Dener ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr); 1590d71dc2bSXiang Huang 160a3c390cfSAlp Dener switch (gn->reg_type) { 161470ec3f8SXiang Huang case BRGN_REGULARIZATION_USER: 162a3c390cfSAlp Dener ierr = (*gn->regularizerhessian)(tao,X,gn->Hreg,gn->reg_hess_ctx);CHKERRQ(ierr); 163a3c390cfSAlp Dener break; 164a1c74439SHansol Suh case BRGN_REGULARIZATION_L2PURE: 165a1c74439SHansol Suh break; 166470ec3f8SXiang Huang case BRGN_REGULARIZATION_L2PROX: 167a3c390cfSAlp Dener break; 168470ec3f8SXiang Huang case BRGN_REGULARIZATION_L1DICT: 1697cea06e1SXiang Huang /* calculate and store diagonal matrix as a vector: diag = epsilon^2 ./ sqrt(x.^2+epsilon^2).^3* --> diag = epsilon^2 ./ sqrt(y.^2+epsilon^2).^3,where y = D*x */ 170a3c390cfSAlp Dener if (gn->D) { 1717cea06e1SXiang Huang ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */ 172a3c390cfSAlp Dener } else { 173a3c390cfSAlp Dener ierr = VecCopy(X,gn->y);CHKERRQ(ierr); 174a3c390cfSAlp Dener } 1757cea06e1SXiang Huang ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr); 1767cea06e1SXiang Huang ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 1777cea06e1SXiang Huang ierr = VecCopy(gn->y_work,gn->diag);CHKERRQ(ierr); /* gn->diag = y.^2+epsilon^2 */ 1787cea06e1SXiang Huang ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr); /* gn->y_work = sqrt(y.^2+epsilon^2) */ 1797cea06e1SXiang Huang ierr = VecPointwiseMult(gn->diag,gn->y_work,gn->diag);CHKERRQ(ierr);/* gn->diag = sqrt(y.^2+epsilon^2).^3 */ 1808ac80d48SXiang Huang ierr = VecReciprocal(gn->diag);CHKERRQ(ierr); 1818ac80d48SXiang Huang ierr = VecScale(gn->diag,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 182a3c390cfSAlp Dener break; 183*cd1c4666STristan Konolige case BRGN_REGULARIZATION_LM: 184*cd1c4666STristan Konolige /* compute diagonal of J^T J */ 185*cd1c4666STristan Konolige ierr = MatGetSize(gn->parent->ls_jac,NULL,&n);CHKERRQ(ierr); 186*cd1c4666STristan Konolige ierr = PetscMalloc1(n,&cnorms);CHKERRQ(ierr); 187*cd1c4666STristan Konolige ierr = MatGetColumnNorms(gn->parent->ls_jac,NORM_2,cnorms);CHKERRQ(ierr); 188*cd1c4666STristan Konolige ierr = MatGetOwnershipRangeColumn(gn->parent->ls_jac,&cstart,&cend);CHKERRQ(ierr); 189*cd1c4666STristan Konolige ierr = VecGetArray(gn->diag,&diag_ary);CHKERRQ(ierr); 190*cd1c4666STristan Konolige for (i = 0; i < cend-cstart; i++) { 191*cd1c4666STristan Konolige diag_ary[i] = cnorms[cstart+i] * cnorms[cstart+i]; 192*cd1c4666STristan Konolige } 193*cd1c4666STristan Konolige ierr = VecRestoreArray(gn->diag,&diag_ary);CHKERRQ(ierr); 194*cd1c4666STristan Konolige ierr = PetscFree(cnorms);CHKERRQ(ierr); 195*cd1c4666STristan Konolige ierr = ComputeDamping(gn);CHKERRQ(ierr); 196*cd1c4666STristan Konolige break; 197a3c390cfSAlp Dener } 198e1e80dc8SAlp Dener PetscFunctionReturn(0); 199e1e80dc8SAlp Dener } 200e1e80dc8SAlp Dener 2018fcddce6SStefano Zampini static PetscErrorCode GNHookFunction(Tao tao,PetscInt iter, void *ctx) 202e1e80dc8SAlp Dener { 2038fcddce6SStefano Zampini TAO_BRGN *gn = (TAO_BRGN *)ctx; 204e1e80dc8SAlp Dener PetscErrorCode ierr; 205e1e80dc8SAlp Dener 206e1e80dc8SAlp Dener PetscFunctionBegin; 207e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */ 208e1e80dc8SAlp Dener gn->parent->nfuncs = tao->nfuncs; 209e1e80dc8SAlp Dener gn->parent->ngrads = tao->ngrads; 210e1e80dc8SAlp Dener gn->parent->nfuncgrads = tao->nfuncgrads; 211e1e80dc8SAlp Dener gn->parent->nhess = tao->nhess; 212e1e80dc8SAlp Dener gn->parent->niter = tao->niter; 213e1e80dc8SAlp Dener gn->parent->ksp_its = tao->ksp_its; 214e1e80dc8SAlp Dener gn->parent->ksp_tot_its = tao->ksp_tot_its; 215*cd1c4666STristan Konolige gn->parent->fc = tao->fc; 216e1e80dc8SAlp Dener ierr = TaoGetConvergedReason(tao,&gn->parent->reason);CHKERRQ(ierr); 217e1e80dc8SAlp Dener /* Update the solution vectors */ 218e1e80dc8SAlp Dener if (iter == 0) { 219e1e80dc8SAlp Dener ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr); 220e1e80dc8SAlp Dener } else { 221e1e80dc8SAlp Dener ierr = VecCopy(tao->solution,gn->x_old);CHKERRQ(ierr); 222e1e80dc8SAlp Dener ierr = VecCopy(tao->solution,gn->parent->solution);CHKERRQ(ierr); 223e1e80dc8SAlp Dener } 224e1e80dc8SAlp Dener /* Update the gradient */ 225e1e80dc8SAlp Dener ierr = VecCopy(tao->gradient,gn->parent->gradient);CHKERRQ(ierr); 226*cd1c4666STristan Konolige 227*cd1c4666STristan Konolige /* Update damping parameter for LM */ 228*cd1c4666STristan Konolige if (gn->reg_type == BRGN_REGULARIZATION_LM) { 229*cd1c4666STristan Konolige if (iter > 0) { 230*cd1c4666STristan Konolige if (gn->fc_old > tao->fc) { 231*cd1c4666STristan Konolige gn->lambda = gn->lambda * gn->downhill_lambda_change; 232*cd1c4666STristan Konolige } else { 233*cd1c4666STristan Konolige /* uphill step */ 234*cd1c4666STristan Konolige gn->lambda = gn->lambda * gn->uphill_lambda_change; 235*cd1c4666STristan Konolige } 236*cd1c4666STristan Konolige } 237*cd1c4666STristan Konolige gn->fc_old = tao->fc; 238*cd1c4666STristan Konolige } 239*cd1c4666STristan Konolige 240e1e80dc8SAlp Dener /* Call general purpose update function */ 241e1e80dc8SAlp Dener if (gn->parent->ops->update) { 2428fcddce6SStefano Zampini ierr = (*gn->parent->ops->update)(gn->parent,gn->parent->niter,gn->parent->user_update);CHKERRQ(ierr); 243737f463aSAlp Dener } 244737f463aSAlp Dener PetscFunctionReturn(0); 245737f463aSAlp Dener } 246737f463aSAlp Dener 247737f463aSAlp Dener static PetscErrorCode TaoSolve_BRGN(Tao tao) 248737f463aSAlp Dener { 249737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 250737f463aSAlp Dener PetscErrorCode ierr; 251737f463aSAlp Dener 252737f463aSAlp Dener PetscFunctionBegin; 253737f463aSAlp Dener ierr = TaoSolve(gn->subsolver);CHKERRQ(ierr); 254e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */ 255e1e80dc8SAlp Dener tao->nfuncs = gn->subsolver->nfuncs; 256e1e80dc8SAlp Dener tao->ngrads = gn->subsolver->ngrads; 257e1e80dc8SAlp Dener tao->nfuncgrads = gn->subsolver->nfuncgrads; 258e1e80dc8SAlp Dener tao->nhess = gn->subsolver->nhess; 259e1e80dc8SAlp Dener tao->niter = gn->subsolver->niter; 260e1e80dc8SAlp Dener tao->ksp_its = gn->subsolver->ksp_its; 261e1e80dc8SAlp Dener tao->ksp_tot_its = gn->subsolver->ksp_tot_its; 262e1e80dc8SAlp Dener ierr = TaoGetConvergedReason(gn->subsolver,&tao->reason);CHKERRQ(ierr); 263e1e80dc8SAlp Dener /* Update vectors */ 264e1e80dc8SAlp Dener ierr = VecCopy(gn->subsolver->solution,tao->solution);CHKERRQ(ierr); 265e1e80dc8SAlp Dener ierr = VecCopy(gn->subsolver->gradient,tao->gradient);CHKERRQ(ierr); 266737f463aSAlp Dener PetscFunctionReturn(0); 267737f463aSAlp Dener } 268737f463aSAlp Dener 269737f463aSAlp Dener static PetscErrorCode TaoSetFromOptions_BRGN(PetscOptionItems *PetscOptionsObject,Tao tao) 270737f463aSAlp Dener { 271737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 272*cd1c4666STristan Konolige TaoLineSearch ls; 273737f463aSAlp Dener PetscErrorCode ierr; 274737f463aSAlp Dener 275737f463aSAlp Dener PetscFunctionBegin; 276a3c390cfSAlp Dener ierr = 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.");CHKERRQ(ierr); 27788fa4459SAlp Dener ierr = PetscOptionsReal("-tao_brgn_regularizer_weight","regularizer weight (default 1e-4)","",gn->lambda,&gn->lambda,NULL);CHKERRQ(ierr); 27888fa4459SAlp Dener ierr = 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);CHKERRQ(ierr); 279*cd1c4666STristan Konolige ierr = 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); 280*cd1c4666STristan Konolige ierr = 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); 281470ec3f8SXiang Huang ierr = PetscOptionsEList("-tao_brgn_regularization_type","regularization type", "",BRGN_REGULARIZATION_TABLE,BRGN_REGULARIZATION_TYPES,BRGN_REGULARIZATION_TABLE[gn->reg_type],&gn->reg_type,NULL);CHKERRQ(ierr); 282737f463aSAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 283*cd1c4666STristan Konolige /* set unit line search direction as the default when using the lm regularizer */ 284*cd1c4666STristan Konolige if (gn->reg_type == BRGN_REGULARIZATION_LM) { 285*cd1c4666STristan Konolige ierr = TaoGetLineSearch(gn->subsolver,&ls);CHKERRQ(ierr); 286*cd1c4666STristan Konolige ierr = TaoLineSearchSetType(ls,TAOLINESEARCHUNIT);CHKERRQ(ierr); 287*cd1c4666STristan Konolige } 288737f463aSAlp Dener ierr = TaoSetFromOptions(gn->subsolver);CHKERRQ(ierr); 289737f463aSAlp Dener PetscFunctionReturn(0); 290737f463aSAlp Dener } 291737f463aSAlp Dener 292737f463aSAlp Dener static PetscErrorCode TaoView_BRGN(Tao tao,PetscViewer viewer) 293737f463aSAlp Dener { 294737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 295737f463aSAlp Dener PetscErrorCode ierr; 296737f463aSAlp Dener 297737f463aSAlp Dener PetscFunctionBegin; 298e1e80dc8SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 299737f463aSAlp Dener ierr = TaoView(gn->subsolver,viewer);CHKERRQ(ierr); 300e1e80dc8SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 301737f463aSAlp Dener PetscFunctionReturn(0); 302737f463aSAlp Dener } 303737f463aSAlp Dener 304737f463aSAlp Dener static PetscErrorCode TaoSetUp_BRGN(Tao tao) 305737f463aSAlp Dener { 306737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 307737f463aSAlp Dener PetscErrorCode ierr; 308737f463aSAlp Dener PetscBool is_bnls,is_bntr,is_bntl; 3098e85b1b3SXiang Huang PetscInt i,n,N,K; /* dict has size K*N*/ 310737f463aSAlp Dener 311737f463aSAlp Dener PetscFunctionBegin; 312737f463aSAlp Dener if (!tao->ls_res) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualRoutine() must be called before setup!"); 313737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNLS,&is_bnls);CHKERRQ(ierr); 314737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTR,&is_bntr);CHKERRQ(ierr); 315737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTL,&is_bntl);CHKERRQ(ierr); 316737f463aSAlp Dener if ((is_bnls || is_bntr || is_bntl) && !tao->ls_jac) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualJacobianRoutine() must be called before setup!"); 317e1e80dc8SAlp Dener if (!tao->gradient) { 318e1e80dc8SAlp Dener ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 319e1e80dc8SAlp Dener } 320737f463aSAlp Dener if (!gn->x_work) { 321737f463aSAlp Dener ierr = VecDuplicate(tao->solution,&gn->x_work);CHKERRQ(ierr); 322737f463aSAlp Dener } 323737f463aSAlp Dener if (!gn->r_work) { 324737f463aSAlp Dener ierr = VecDuplicate(tao->ls_res,&gn->r_work);CHKERRQ(ierr); 325737f463aSAlp Dener } 326e1e80dc8SAlp Dener if (!gn->x_old) { 327e1e80dc8SAlp Dener ierr = VecDuplicate(tao->solution,&gn->x_old);CHKERRQ(ierr); 328e1e80dc8SAlp Dener ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr); 329e1e80dc8SAlp Dener } 3307cea06e1SXiang Huang 331470ec3f8SXiang Huang if (BRGN_REGULARIZATION_L1DICT == gn->reg_type) { 33230eeff36SXiang Huang if (gn->D) { 33330eeff36SXiang Huang ierr = MatGetSize(gn->D,&K,&N);CHKERRQ(ierr); /* Shell matrices still must have sizes defined. K = N for identity matrix, K=N-1 or N for gradient matrix */ 33430eeff36SXiang Huang } else { 33530eeff36SXiang Huang ierr = VecGetSize(tao->solution,&K);CHKERRQ(ierr); /* If user does not setup dict matrix, use identiy matrix, K=N */ 33630eeff36SXiang Huang } 3377cea06e1SXiang Huang if (!gn->y) { 3387cea06e1SXiang Huang ierr = VecCreate(PETSC_COMM_SELF,&gn->y);CHKERRQ(ierr); 3398e85b1b3SXiang Huang ierr = VecSetSizes(gn->y,PETSC_DECIDE,K);CHKERRQ(ierr); 3407cea06e1SXiang Huang ierr = VecSetFromOptions(gn->y);CHKERRQ(ierr); 3417cea06e1SXiang Huang ierr = VecSet(gn->y,0.0);CHKERRQ(ierr); 3427cea06e1SXiang Huang 3437cea06e1SXiang Huang } 3447cea06e1SXiang Huang if (!gn->y_work) { 3457cea06e1SXiang Huang ierr = VecDuplicate(gn->y,&gn->y_work);CHKERRQ(ierr); 3467cea06e1SXiang Huang } 3478ac80d48SXiang Huang if (!gn->diag) { 3487cea06e1SXiang Huang ierr = VecDuplicate(gn->y,&gn->diag);CHKERRQ(ierr); 3498ac80d48SXiang Huang ierr = VecSet(gn->diag,0.0);CHKERRQ(ierr); 3508ac80d48SXiang Huang } 35130eeff36SXiang Huang } 352*cd1c4666STristan Konolige if (BRGN_REGULARIZATION_LM == gn->reg_type) { 353*cd1c4666STristan Konolige if (!gn->diag) { 354*cd1c4666STristan Konolige ierr = MatCreateVecs(gn->parent->ls_jac,&gn->diag,NULL);CHKERRQ(ierr); 355*cd1c4666STristan Konolige } 356*cd1c4666STristan Konolige if (!gn->damping) { 357*cd1c4666STristan Konolige ierr = MatCreateVecs(gn->parent->ls_jac,&gn->damping,NULL);CHKERRQ(ierr); 358*cd1c4666STristan Konolige } 359*cd1c4666STristan Konolige } 3600d71dc2bSXiang Huang 361e1e80dc8SAlp Dener if (!tao->setupcalled) { 362737f463aSAlp Dener /* Hessian setup */ 3638e85b1b3SXiang Huang ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 3648e85b1b3SXiang Huang ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 3658e85b1b3SXiang Huang ierr = MatSetSizes(gn->H,n,n,N,N);CHKERRQ(ierr); 366737f463aSAlp Dener ierr = MatSetType(gn->H,MATSHELL);CHKERRQ(ierr); 367737f463aSAlp Dener ierr = MatSetUp(gn->H);CHKERRQ(ierr); 368737f463aSAlp Dener ierr = MatShellSetOperation(gn->H,MATOP_MULT,(void (*)(void))GNHessianProd);CHKERRQ(ierr); 369737f463aSAlp Dener ierr = MatShellSetContext(gn->H,(void*)gn);CHKERRQ(ierr); 3708e85b1b3SXiang Huang /* Subsolver setup,include initial vector and dicttionary D */ 371e1e80dc8SAlp Dener ierr = TaoSetUpdate(gn->subsolver,GNHookFunction,(void*)gn);CHKERRQ(ierr); 372737f463aSAlp Dener ierr = TaoSetInitialVector(gn->subsolver,tao->solution);CHKERRQ(ierr); 373737f463aSAlp Dener if (tao->bounded) { 374737f463aSAlp Dener ierr = TaoSetVariableBounds(gn->subsolver,tao->XL,tao->XU);CHKERRQ(ierr); 375737f463aSAlp Dener } 376737f463aSAlp Dener ierr = TaoSetResidualRoutine(gn->subsolver,tao->ls_res,tao->ops->computeresidual,tao->user_lsresP);CHKERRQ(ierr); 3774ffbe8acSAlp Dener ierr = TaoSetJacobianResidualRoutine(gn->subsolver,tao->ls_jac,tao->ls_jac,tao->ops->computeresidualjacobian,tao->user_lsjacP);CHKERRQ(ierr); 378737f463aSAlp Dener ierr = TaoSetObjectiveAndGradientRoutine(gn->subsolver,GNObjectiveGradientEval,(void*)gn);CHKERRQ(ierr); 379737f463aSAlp Dener ierr = TaoSetHessianRoutine(gn->subsolver,gn->H,gn->H,GNComputeHessian,(void*)gn);CHKERRQ(ierr); 380e1e80dc8SAlp Dener /* Propagate some options down */ 381e1e80dc8SAlp Dener ierr = TaoSetTolerances(gn->subsolver,tao->gatol,tao->grtol,tao->gttol);CHKERRQ(ierr); 382e1e80dc8SAlp Dener ierr = TaoSetMaximumIterations(gn->subsolver,tao->max_it);CHKERRQ(ierr); 383e1e80dc8SAlp Dener ierr = TaoSetMaximumFunctionEvaluations(gn->subsolver,tao->max_funcs);CHKERRQ(ierr); 384737f463aSAlp Dener for (i=0; i<tao->numbermonitors; ++i) { 385737f463aSAlp Dener ierr = TaoSetMonitor(gn->subsolver,tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i]);CHKERRQ(ierr); 386737f463aSAlp Dener ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr); 387737f463aSAlp Dener } 388737f463aSAlp Dener ierr = TaoSetUp(gn->subsolver);CHKERRQ(ierr); 389e1e80dc8SAlp Dener } 390737f463aSAlp Dener PetscFunctionReturn(0); 391737f463aSAlp Dener } 392737f463aSAlp Dener 393737f463aSAlp Dener static PetscErrorCode TaoDestroy_BRGN(Tao tao) 394737f463aSAlp Dener { 395737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 396737f463aSAlp Dener PetscErrorCode ierr; 397737f463aSAlp Dener 398737f463aSAlp Dener PetscFunctionBegin; 399737f463aSAlp Dener if (tao->setupcalled) { 400e1e80dc8SAlp Dener ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr); 401737f463aSAlp Dener ierr = VecDestroy(&gn->x_work);CHKERRQ(ierr); 402737f463aSAlp Dener ierr = VecDestroy(&gn->r_work);CHKERRQ(ierr); 403e1e80dc8SAlp Dener ierr = VecDestroy(&gn->x_old);CHKERRQ(ierr); 4048ac80d48SXiang Huang ierr = VecDestroy(&gn->diag);CHKERRQ(ierr); 4057cea06e1SXiang Huang ierr = VecDestroy(&gn->y);CHKERRQ(ierr); 4067cea06e1SXiang Huang ierr = VecDestroy(&gn->y_work);CHKERRQ(ierr); 407737f463aSAlp Dener } 408*cd1c4666STristan Konolige ierr = VecDestroy(&gn->damping);CHKERRQ(ierr); 409*cd1c4666STristan Konolige ierr = VecDestroy(&gn->diag);CHKERRQ(ierr); 410737f463aSAlp Dener ierr = MatDestroy(&gn->H);CHKERRQ(ierr); 4117cea06e1SXiang Huang ierr = MatDestroy(&gn->D);CHKERRQ(ierr); 412463fc0ecSAlp Dener ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr); 413737f463aSAlp Dener ierr = TaoDestroy(&gn->subsolver);CHKERRQ(ierr); 414e1e80dc8SAlp Dener gn->parent = NULL; 415737f463aSAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 416737f463aSAlp Dener PetscFunctionReturn(0); 417737f463aSAlp Dener } 418737f463aSAlp Dener 4193850be85SAlp Dener /*MC 4203850be85SAlp Dener TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares 4213850be85SAlp Dener problems with bound constraints. This algorithm is a thin wrapper around TAOBNTL 422463fc0ecSAlp Dener that constructs the Gauss-Newton problem with the user-provided least-squares 42360bb7533SHansol Suh residual and Jacobian. The algorithm offers an L2-norm ("l2pure"), L2-norm proximal point ("l2prox") 42460bb7533SHansol Suh regularizer, and L1-norm dictionary regularizer ("l1dict"), where we approximate the 42501b716f5SXiang Huang L1-norm ||x||_1 by sum_i(sqrt(x_i^2+epsilon^2)-epsilon) with a small positive number epsilon. 426*cd1c4666STristan Konolige Also offered is the "lm" regularizer which uses a scaled diagonal of J^T J. 427*cd1c4666STristan Konolige With the "lm" regularizer, BRGN is a Levenberg-Marquardt optimizer. 42801b716f5SXiang Huang The user can also provide own regularization function. 4293850be85SAlp Dener 4303850be85SAlp Dener Options Database Keys: 431*cd1c4666STristan Konolige + -tao_brgn_regularization_type - regularization type ("user", "l2prox", "l2pure", "l1dict", "lm") (default "l2prox") 432c061e8e2SXiang Huang . -tao_brgn_regularizer_weight - regularizer weight (default 1e-4) 433c061e8e2SXiang Huang - -tao_brgn_l1_smooth_epsilon - L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6) 4343850be85SAlp Dener 4353850be85SAlp Dener Level: beginner 4363850be85SAlp Dener M*/ 437737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao) 438737f463aSAlp Dener { 439737f463aSAlp Dener TAO_BRGN *gn; 440737f463aSAlp Dener PetscErrorCode ierr; 441737f463aSAlp Dener 442737f463aSAlp Dener PetscFunctionBegin; 443737f463aSAlp Dener ierr = PetscNewLog(tao,&gn);CHKERRQ(ierr); 444737f463aSAlp Dener 445737f463aSAlp Dener tao->ops->destroy = TaoDestroy_BRGN; 446737f463aSAlp Dener tao->ops->setup = TaoSetUp_BRGN; 447737f463aSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BRGN; 448737f463aSAlp Dener tao->ops->view = TaoView_BRGN; 449737f463aSAlp Dener tao->ops->solve = TaoSolve_BRGN; 450737f463aSAlp Dener 451737f463aSAlp Dener tao->data = (void*)gn; 452d8bf7057SXiang Huang gn->reg_type = BRGN_REGULARIZATION_L2PROX; 453e1e80dc8SAlp Dener gn->lambda = 1e-4; 4548ac80d48SXiang Huang gn->epsilon = 1e-6; 455*cd1c4666STristan Konolige gn->downhill_lambda_change = 1./5.; 456*cd1c4666STristan Konolige gn->uphill_lambda_change = 1.5; 457e1e80dc8SAlp Dener gn->parent = tao; 458737f463aSAlp Dener 459737f463aSAlp Dener ierr = MatCreate(PetscObjectComm((PetscObject)tao),&gn->H);CHKERRQ(ierr); 460737f463aSAlp Dener ierr = MatSetOptionsPrefix(gn->H,"tao_brgn_hessian_");CHKERRQ(ierr); 461737f463aSAlp Dener 462737f463aSAlp Dener ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&gn->subsolver);CHKERRQ(ierr); 463737f463aSAlp Dener ierr = TaoSetType(gn->subsolver,TAOBNLS);CHKERRQ(ierr); 464737f463aSAlp Dener ierr = TaoSetOptionsPrefix(gn->subsolver,"tao_brgn_subsolver_");CHKERRQ(ierr); 465e1e80dc8SAlp Dener PetscFunctionReturn(0); 466e1e80dc8SAlp Dener } 467e1e80dc8SAlp Dener 468463fc0ecSAlp Dener /*@ 469e1e80dc8SAlp Dener TaoBRGNGetSubsolver - Get the pointer to the subsolver inside BRGN 470e1e80dc8SAlp Dener 471e1e80dc8SAlp Dener Collective on Tao 472e1e80dc8SAlp Dener 473463fc0ecSAlp Dener Level: advanced 474e1e80dc8SAlp Dener 475e1e80dc8SAlp Dener Input Parameters: 476e1e80dc8SAlp Dener + tao - the Tao solver context 477e1e80dc8SAlp Dener - subsolver - the Tao sub-solver context 478e1e80dc8SAlp Dener @*/ 479e1e80dc8SAlp Dener PetscErrorCode TaoBRGNGetSubsolver(Tao tao,Tao *subsolver) 480e1e80dc8SAlp Dener { 481e1e80dc8SAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 482e1e80dc8SAlp Dener 483e1e80dc8SAlp Dener PetscFunctionBegin; 484e1e80dc8SAlp Dener *subsolver = gn->subsolver; 485737f463aSAlp Dener PetscFunctionReturn(0); 486737f463aSAlp Dener } 487737f463aSAlp Dener 488463fc0ecSAlp Dener /*@ 489463fc0ecSAlp Dener TaoBRGNSetRegularizerWeight - Set the regularizer weight for the Gauss-Newton least-squares algorithm 490737f463aSAlp Dener 491737f463aSAlp Dener Collective on Tao 492737f463aSAlp Dener 493737f463aSAlp Dener Input Parameters: 494737f463aSAlp Dener + tao - the Tao solver context 4958e85b1b3SXiang Huang - lambda - L1-norm regularizer weight 496463fc0ecSAlp Dener 497463fc0ecSAlp Dener Level: beginner 498737f463aSAlp Dener @*/ 499a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao,PetscReal lambda) 500737f463aSAlp Dener { 501737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 502737f463aSAlp Dener 5038ac80d48SXiang Huang /* Initialize lambda here */ 5040d71dc2bSXiang Huang 505737f463aSAlp Dener PetscFunctionBegin; 506737f463aSAlp Dener gn->lambda = lambda; 507737f463aSAlp Dener PetscFunctionReturn(0); 508737f463aSAlp Dener } 5090d71dc2bSXiang Huang 510463fc0ecSAlp Dener /*@ 5118ac80d48SXiang Huang TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm 5128ac80d48SXiang Huang 5138ac80d48SXiang Huang Collective on Tao 5148ac80d48SXiang Huang 5158ac80d48SXiang Huang Input Parameters: 5168ac80d48SXiang Huang + tao - the Tao solver context 5178ac80d48SXiang Huang - epsilon - L1-norm smooth approximation parameter 518463fc0ecSAlp Dener 519463fc0ecSAlp Dener Level: advanced 5208ac80d48SXiang Huang @*/ 5218ac80d48SXiang Huang PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao,PetscReal epsilon) 5228ac80d48SXiang Huang { 5238ac80d48SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)tao->data; 5248ac80d48SXiang Huang 5258ac80d48SXiang Huang /* Initialize epsilon here */ 5268ac80d48SXiang Huang 5278ac80d48SXiang Huang PetscFunctionBegin; 5288ac80d48SXiang Huang gn->epsilon = epsilon; 5298ac80d48SXiang Huang PetscFunctionReturn(0); 5308ac80d48SXiang Huang } 5318e85b1b3SXiang Huang 532463fc0ecSAlp Dener /*@ 5338e85b1b3SXiang Huang TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem) 5348e85b1b3SXiang Huang 5358e85b1b3SXiang Huang Input Parameters: 5368e85b1b3SXiang Huang + tao - the Tao context 537a2b725a8SWilliam Gropp - dict - the user specified dictionary matrix. We allow to set a null dictionary, which means identity matrix by default 5388e85b1b3SXiang Huang 539463fc0ecSAlp Dener Level: advanced 5408e85b1b3SXiang Huang @*/ 5418e85b1b3SXiang Huang PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao,Mat dict) 5428e85b1b3SXiang Huang { 5438e85b1b3SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)tao->data; 5448e85b1b3SXiang Huang PetscErrorCode ierr; 5458e85b1b3SXiang Huang PetscFunctionBegin; 5468e85b1b3SXiang Huang PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 5478e85b1b3SXiang Huang if (dict) { 5488e85b1b3SXiang Huang PetscValidHeaderSpecific(dict,MAT_CLASSID,2); 549a3c390cfSAlp Dener PetscCheckSameComm(tao,1,dict,2); 5508e85b1b3SXiang Huang ierr = PetscObjectReference((PetscObject)dict);CHKERRQ(ierr); 5518e85b1b3SXiang Huang } 5528e85b1b3SXiang Huang ierr = MatDestroy(&gn->D);CHKERRQ(ierr); 5531fc140a9SXiang Huang gn->D = dict; 5548e85b1b3SXiang Huang PetscFunctionReturn(0); 5558e85b1b3SXiang Huang } 5568e85b1b3SXiang Huang 557a3c390cfSAlp Dener /*@C 558463fc0ecSAlp Dener TaoBRGNSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back 559463fc0ecSAlp Dener function into the algorithm. 560463fc0ecSAlp Dener 561463fc0ecSAlp Dener Input Parameters: 562463fc0ecSAlp Dener + tao - the Tao context 563463fc0ecSAlp Dener . func - function pointer for the regularizer value and gradient evaluation 564463fc0ecSAlp Dener - ctx - user context for the regularizer 565463fc0ecSAlp Dener 566463fc0ecSAlp Dener Level: advanced 567a3c390cfSAlp Dener @*/ 568a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (*func)(Tao,Vec,PetscReal *,Vec,void*),void *ctx) 569a3c390cfSAlp Dener { 570a3c390cfSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 571a3c390cfSAlp Dener 572a3c390cfSAlp Dener PetscFunctionBegin; 573a3c390cfSAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 574a3c390cfSAlp Dener if (ctx) { 575a3c390cfSAlp Dener gn->reg_obj_ctx = ctx; 576a3c390cfSAlp Dener } 577a3c390cfSAlp Dener if (func) { 578a3c390cfSAlp Dener gn->regularizerobjandgrad = func; 579a3c390cfSAlp Dener } 580a3c390cfSAlp Dener PetscFunctionReturn(0); 581a3c390cfSAlp Dener } 582a3c390cfSAlp Dener 583a3c390cfSAlp Dener /*@C 584463fc0ecSAlp Dener TaoBRGNSetRegularizerHessianRoutine - Sets the user-defined regularizer call-back 585463fc0ecSAlp Dener function into the algorithm. 586463fc0ecSAlp Dener 587463fc0ecSAlp Dener Input Parameters: 588463fc0ecSAlp Dener + tao - the Tao context 589463fc0ecSAlp Dener . Hreg - user-created matrix for the Hessian of the regularization term 590463fc0ecSAlp Dener . func - function pointer for the regularizer Hessian evaluation 591463fc0ecSAlp Dener - ctx - user context for the regularizer Hessian 592463fc0ecSAlp Dener 593463fc0ecSAlp Dener Level: advanced 594a3c390cfSAlp Dener @*/ 595a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao,Mat Hreg,PetscErrorCode (*func)(Tao,Vec,Mat,void*),void *ctx) 596a3c390cfSAlp Dener { 597a3c390cfSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 598a3c390cfSAlp Dener PetscErrorCode ierr; 599a3c390cfSAlp Dener 600a3c390cfSAlp Dener PetscFunctionBegin; 601a3c390cfSAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 602a3c390cfSAlp Dener if (Hreg) { 603a3c390cfSAlp Dener PetscValidHeaderSpecific(Hreg,MAT_CLASSID,2); 604a3c390cfSAlp Dener PetscCheckSameComm(tao,1,Hreg,2); 60501b716f5SXiang Huang } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONG,"NULL Hessian detected! User must provide valid Hessian for the regularizer."); 606a3c390cfSAlp Dener if (ctx) { 607a3c390cfSAlp Dener gn->reg_hess_ctx = ctx; 608a3c390cfSAlp Dener } 609a3c390cfSAlp Dener if (func) { 610a3c390cfSAlp Dener gn->regularizerhessian = func; 611a3c390cfSAlp Dener } 612a3c390cfSAlp Dener if (Hreg) { 613a3c390cfSAlp Dener ierr = PetscObjectReference((PetscObject)Hreg);CHKERRQ(ierr); 614a3c390cfSAlp Dener ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr); 615a3c390cfSAlp Dener gn->Hreg = Hreg; 616a3c390cfSAlp Dener } 617a3c390cfSAlp Dener PetscFunctionReturn(0); 618a3c390cfSAlp Dener } 619