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 5470ec3f8SXiang Huang #define BRGN_REGULARIZATION_L1DICT 2 6470ec3f8SXiang Huang #define BRGN_REGULARIZATION_TYPES 3 7a3c390cfSAlp Dener 8470ec3f8SXiang Huang static const char *BRGN_REGULARIZATION_TABLE[64] = {"user","l2prox","l1dict"}; 9a3c390cfSAlp Dener 100d71dc2bSXiang Huang static PetscErrorCode GNHessianProd(Mat H,Vec in,Vec out) 110d71dc2bSXiang Huang { 120d71dc2bSXiang Huang TAO_BRGN *gn; 130d71dc2bSXiang Huang PetscErrorCode ierr; 140d71dc2bSXiang Huang 150d71dc2bSXiang Huang PetscFunctionBegin; 160d71dc2bSXiang Huang ierr = MatShellGetContext(H,&gn);CHKERRQ(ierr); 170d71dc2bSXiang Huang ierr = MatMult(gn->subsolver->ls_jac,in,gn->r_work);CHKERRQ(ierr); 180d71dc2bSXiang Huang ierr = MatMultTranspose(gn->subsolver->ls_jac,gn->r_work,out);CHKERRQ(ierr); 19a3c390cfSAlp Dener switch (gn->reg_type) { 20470ec3f8SXiang Huang case BRGN_REGULARIZATION_USER: 21a3c390cfSAlp Dener ierr = MatMult(gn->Hreg,in,gn->x_work);CHKERRQ(ierr); 228ac80d48SXiang Huang ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr); 23a3c390cfSAlp Dener break; 24470ec3f8SXiang Huang case BRGN_REGULARIZATION_L2PROX: 25a3c390cfSAlp Dener ierr = VecAXPY(out,gn->lambda,in);CHKERRQ(ierr); 26a3c390cfSAlp Dener break; 27470ec3f8SXiang Huang case BRGN_REGULARIZATION_L1DICT: 28a3c390cfSAlp Dener /* out = out + lambda*D'*(diag.*(D*in)) */ 29a3c390cfSAlp Dener if (gn->D) { 30a3c390cfSAlp Dener ierr = MatMult(gn->D,in,gn->y);CHKERRQ(ierr);/* y = D*in */ 31a3c390cfSAlp Dener } else { 32a3c390cfSAlp Dener ierr = VecCopy(in,gn->y);CHKERRQ(ierr); 33a3c390cfSAlp Dener } 34a3c390cfSAlp 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 */ 35a3c390cfSAlp Dener if (gn->D) { 36a3c390cfSAlp Dener ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr); /* x_work = D'*(diag.*(D*in)) */ 37a3c390cfSAlp Dener } else { 38a3c390cfSAlp Dener ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr); 39a3c390cfSAlp Dener } 40a3c390cfSAlp Dener ierr = VecAXPY(out,gn->lambda,gn->x_work);CHKERRQ(ierr); 41a3c390cfSAlp Dener break; 42a3c390cfSAlp Dener } 430d71dc2bSXiang Huang PetscFunctionReturn(0); 440d71dc2bSXiang Huang } 450d71dc2bSXiang Huang 460d71dc2bSXiang Huang static PetscErrorCode GNObjectiveGradientEval(Tao tao,Vec X,PetscReal *fcn,Vec G,void *ptr) 470d71dc2bSXiang Huang { 480d71dc2bSXiang Huang TAO_BRGN *gn = (TAO_BRGN *)ptr; 498e85b1b3SXiang Huang PetscInt K; /* dimension of D*X */ 507cea06e1SXiang Huang PetscScalar yESum; 510d71dc2bSXiang Huang PetscErrorCode ierr; 52a3c390cfSAlp Dener PetscReal f_reg; 530d71dc2bSXiang Huang 540d71dc2bSXiang Huang PetscFunctionBegin; 558e85b1b3SXiang Huang /* compute objective *fcn*/ 56a3c390cfSAlp Dener /* compute first term 0.5*||ls_res||_2^2 */ 570d71dc2bSXiang Huang ierr = TaoComputeResidual(tao,X,tao->ls_res);CHKERRQ(ierr); 58a3c390cfSAlp Dener ierr = VecDot(tao->ls_res,tao->ls_res,fcn);CHKERRQ(ierr); 59a3c390cfSAlp Dener *fcn *= 0.5; 60a3c390cfSAlp Dener /* compute gradient G */ 61a3c390cfSAlp Dener ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr); 62a3c390cfSAlp Dener ierr = MatMultTranspose(tao->ls_jac,tao->ls_res,G);CHKERRQ(ierr); 63a3c390cfSAlp Dener /* add the regularization contribution */ 64a3c390cfSAlp Dener switch (gn->reg_type) { 65470ec3f8SXiang Huang case BRGN_REGULARIZATION_USER: 66a3c390cfSAlp Dener ierr = (*gn->regularizerobjandgrad)(tao,X,&f_reg,gn->x_work,gn->reg_obj_ctx);CHKERRQ(ierr); 67a3c390cfSAlp Dener *fcn += gn->lambda*f_reg; 68a3c390cfSAlp Dener ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr); 69a3c390cfSAlp Dener break; 70470ec3f8SXiang Huang case BRGN_REGULARIZATION_L2PROX: 711fc140a9SXiang Huang /* compute f = f + lambda*0.5*(xk - xkm1)'*(xk - xkm1) */ 721fc140a9SXiang Huang ierr = VecAXPBYPCZ(gn->x_work,1.0,-1.0,0.0,X,gn->x_old);CHKERRQ(ierr); 73a3c390cfSAlp Dener ierr = VecDot(gn->x_work,gn->x_work,&f_reg);CHKERRQ(ierr); 74a3c390cfSAlp Dener *fcn += gn->lambda*0.5*f_reg; 75a3c390cfSAlp Dener /* compute G = G + lambda*(xk - xkm1) */ 76a3c390cfSAlp Dener ierr = VecAXPBYPCZ(G,gn->lambda,-gn->lambda,1.0,X,gn->x_old);CHKERRQ(ierr); 77a3c390cfSAlp Dener break; 78470ec3f8SXiang Huang case BRGN_REGULARIZATION_L1DICT: 79a3c390cfSAlp Dener /* compute f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x*/ 80a3c390cfSAlp Dener if (gn->D) { 817cea06e1SXiang Huang ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */ 82a3c390cfSAlp Dener } else { 83a3c390cfSAlp Dener ierr = VecCopy(X,gn->y);CHKERRQ(ierr); 84a3c390cfSAlp Dener } 857cea06e1SXiang Huang ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr); 867cea06e1SXiang Huang ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 877cea06e1SXiang Huang ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr); /* gn->y_work = sqrt(y.^2+epsilon^2) */ 887cea06e1SXiang Huang ierr = VecSum(gn->y_work,&yESum);CHKERRQ(ierr);CHKERRQ(ierr); 898e85b1b3SXiang Huang ierr = VecGetSize(gn->y,&K);CHKERRQ(ierr); 90a3c390cfSAlp Dener *fcn += gn->lambda*(yESum - K*gn->epsilon); 917cea06e1SXiang Huang /* compute G = G + lambda*D'*(y./sqrt(y.^2+epsilon^2)),where y = D*x */ 927cea06e1SXiang Huang ierr = VecPointwiseDivide(gn->y_work,gn->y,gn->y_work);CHKERRQ(ierr); /* reuse y_work = y./sqrt(y.^2+epsilon^2) */ 93a3c390cfSAlp Dener if (gn->D) { 947cea06e1SXiang Huang ierr = MatMultTranspose(gn->D,gn->y_work,gn->x_work);CHKERRQ(ierr); 95a3c390cfSAlp Dener } else { 96a3c390cfSAlp Dener ierr = VecCopy(gn->y_work,gn->x_work);CHKERRQ(ierr); 97a3c390cfSAlp Dener } 988ac80d48SXiang Huang ierr = VecAXPY(G,gn->lambda,gn->x_work);CHKERRQ(ierr); 99a3c390cfSAlp Dener break; 100a3c390cfSAlp Dener } 1010d71dc2bSXiang Huang PetscFunctionReturn(0); 1020d71dc2bSXiang Huang } 1030d71dc2bSXiang Huang 104737f463aSAlp Dener static PetscErrorCode GNComputeHessian(Tao tao,Vec X,Mat H,Mat Hpre,void *ptr) 105737f463aSAlp Dener { 1068ac80d48SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)ptr; 107737f463aSAlp Dener PetscErrorCode ierr; 108737f463aSAlp Dener 109737f463aSAlp Dener PetscFunctionBegin; 110e1e80dc8SAlp Dener ierr = TaoComputeResidualJacobian(tao,X,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr); 1110d71dc2bSXiang Huang 112a3c390cfSAlp Dener switch (gn->reg_type) { 113470ec3f8SXiang Huang case BRGN_REGULARIZATION_USER: 114a3c390cfSAlp Dener ierr = (*gn->regularizerhessian)(tao,X,gn->Hreg,gn->reg_hess_ctx);CHKERRQ(ierr); 115a3c390cfSAlp Dener break; 116470ec3f8SXiang Huang case BRGN_REGULARIZATION_L2PROX: 117a3c390cfSAlp Dener break; 118470ec3f8SXiang Huang case BRGN_REGULARIZATION_L1DICT: 1197cea06e1SXiang 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 */ 120a3c390cfSAlp Dener if (gn->D) { 1217cea06e1SXiang Huang ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */ 122a3c390cfSAlp Dener } else { 123a3c390cfSAlp Dener ierr = VecCopy(X,gn->y);CHKERRQ(ierr); 124a3c390cfSAlp Dener } 1257cea06e1SXiang Huang ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr); 1267cea06e1SXiang Huang ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 1277cea06e1SXiang Huang ierr = VecCopy(gn->y_work,gn->diag);CHKERRQ(ierr); /* gn->diag = y.^2+epsilon^2 */ 1287cea06e1SXiang Huang ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr); /* gn->y_work = sqrt(y.^2+epsilon^2) */ 1297cea06e1SXiang Huang ierr = VecPointwiseMult(gn->diag,gn->y_work,gn->diag);CHKERRQ(ierr);/* gn->diag = sqrt(y.^2+epsilon^2).^3 */ 1308ac80d48SXiang Huang ierr = VecReciprocal(gn->diag);CHKERRQ(ierr); 1318ac80d48SXiang Huang ierr = VecScale(gn->diag,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 132a3c390cfSAlp Dener break; 133a3c390cfSAlp Dener } 134e1e80dc8SAlp Dener PetscFunctionReturn(0); 135e1e80dc8SAlp Dener } 136e1e80dc8SAlp Dener 137*8fcddce6SStefano Zampini static PetscErrorCode GNHookFunction(Tao tao,PetscInt iter, void *ctx) 138e1e80dc8SAlp Dener { 139*8fcddce6SStefano Zampini TAO_BRGN *gn = (TAO_BRGN *)ctx; 140e1e80dc8SAlp Dener PetscErrorCode ierr; 141e1e80dc8SAlp Dener 142e1e80dc8SAlp Dener PetscFunctionBegin; 143e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */ 144e1e80dc8SAlp Dener gn->parent->nfuncs = tao->nfuncs; 145e1e80dc8SAlp Dener gn->parent->ngrads = tao->ngrads; 146e1e80dc8SAlp Dener gn->parent->nfuncgrads = tao->nfuncgrads; 147e1e80dc8SAlp Dener gn->parent->nhess = tao->nhess; 148e1e80dc8SAlp Dener gn->parent->niter = tao->niter; 149e1e80dc8SAlp Dener gn->parent->ksp_its = tao->ksp_its; 150e1e80dc8SAlp Dener gn->parent->ksp_tot_its = tao->ksp_tot_its; 151e1e80dc8SAlp Dener ierr = TaoGetConvergedReason(tao,&gn->parent->reason);CHKERRQ(ierr); 152e1e80dc8SAlp Dener /* Update the solution vectors */ 153e1e80dc8SAlp Dener if (iter == 0) { 154e1e80dc8SAlp Dener ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr); 155e1e80dc8SAlp Dener } else { 156e1e80dc8SAlp Dener ierr = VecCopy(tao->solution,gn->x_old);CHKERRQ(ierr); 157e1e80dc8SAlp Dener ierr = VecCopy(tao->solution,gn->parent->solution);CHKERRQ(ierr); 158e1e80dc8SAlp Dener } 159e1e80dc8SAlp Dener /* Update the gradient */ 160e1e80dc8SAlp Dener ierr = VecCopy(tao->gradient,gn->parent->gradient);CHKERRQ(ierr); 161e1e80dc8SAlp Dener /* Call general purpose update function */ 162e1e80dc8SAlp Dener if (gn->parent->ops->update) { 163*8fcddce6SStefano Zampini ierr = (*gn->parent->ops->update)(gn->parent,gn->parent->niter,gn->parent->user_update);CHKERRQ(ierr); 164737f463aSAlp Dener } 165737f463aSAlp Dener PetscFunctionReturn(0); 166737f463aSAlp Dener } 167737f463aSAlp Dener 168737f463aSAlp Dener static PetscErrorCode TaoSolve_BRGN(Tao tao) 169737f463aSAlp Dener { 170737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 171737f463aSAlp Dener PetscErrorCode ierr; 172737f463aSAlp Dener 173737f463aSAlp Dener PetscFunctionBegin; 174737f463aSAlp Dener ierr = TaoSolve(gn->subsolver);CHKERRQ(ierr); 175e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */ 176e1e80dc8SAlp Dener tao->nfuncs = gn->subsolver->nfuncs; 177e1e80dc8SAlp Dener tao->ngrads = gn->subsolver->ngrads; 178e1e80dc8SAlp Dener tao->nfuncgrads = gn->subsolver->nfuncgrads; 179e1e80dc8SAlp Dener tao->nhess = gn->subsolver->nhess; 180e1e80dc8SAlp Dener tao->niter = gn->subsolver->niter; 181e1e80dc8SAlp Dener tao->ksp_its = gn->subsolver->ksp_its; 182e1e80dc8SAlp Dener tao->ksp_tot_its = gn->subsolver->ksp_tot_its; 183e1e80dc8SAlp Dener ierr = TaoGetConvergedReason(gn->subsolver,&tao->reason);CHKERRQ(ierr); 184e1e80dc8SAlp Dener /* Update vectors */ 185e1e80dc8SAlp Dener ierr = VecCopy(gn->subsolver->solution,tao->solution);CHKERRQ(ierr); 186e1e80dc8SAlp Dener ierr = VecCopy(gn->subsolver->gradient,tao->gradient);CHKERRQ(ierr); 187737f463aSAlp Dener PetscFunctionReturn(0); 188737f463aSAlp Dener } 189737f463aSAlp Dener 190737f463aSAlp Dener static PetscErrorCode TaoSetFromOptions_BRGN(PetscOptionItems *PetscOptionsObject,Tao tao) 191737f463aSAlp Dener { 192737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 193737f463aSAlp Dener PetscErrorCode ierr; 194737f463aSAlp Dener 195737f463aSAlp Dener PetscFunctionBegin; 196a3c390cfSAlp 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); 19788fa4459SAlp Dener ierr = PetscOptionsReal("-tao_brgn_regularizer_weight","regularizer weight (default 1e-4)","",gn->lambda,&gn->lambda,NULL);CHKERRQ(ierr); 19888fa4459SAlp 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); 199470ec3f8SXiang 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); 200737f463aSAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 201737f463aSAlp Dener ierr = TaoSetFromOptions(gn->subsolver);CHKERRQ(ierr); 202737f463aSAlp Dener PetscFunctionReturn(0); 203737f463aSAlp Dener } 204737f463aSAlp Dener 205737f463aSAlp Dener static PetscErrorCode TaoView_BRGN(Tao tao,PetscViewer viewer) 206737f463aSAlp Dener { 207737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 208737f463aSAlp Dener PetscErrorCode ierr; 209737f463aSAlp Dener 210737f463aSAlp Dener PetscFunctionBegin; 211e1e80dc8SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 212737f463aSAlp Dener ierr = TaoView(gn->subsolver,viewer);CHKERRQ(ierr); 213e1e80dc8SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 214737f463aSAlp Dener PetscFunctionReturn(0); 215737f463aSAlp Dener } 216737f463aSAlp Dener 217737f463aSAlp Dener static PetscErrorCode TaoSetUp_BRGN(Tao tao) 218737f463aSAlp Dener { 219737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 220737f463aSAlp Dener PetscErrorCode ierr; 221737f463aSAlp Dener PetscBool is_bnls,is_bntr,is_bntl; 2228e85b1b3SXiang Huang PetscInt i,n,N,K; /* dict has size K*N*/ 223737f463aSAlp Dener 224737f463aSAlp Dener PetscFunctionBegin; 225737f463aSAlp Dener if (!tao->ls_res) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualRoutine() must be called before setup!"); 226737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNLS,&is_bnls);CHKERRQ(ierr); 227737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTR,&is_bntr);CHKERRQ(ierr); 228737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTL,&is_bntl);CHKERRQ(ierr); 229737f463aSAlp Dener if ((is_bnls || is_bntr || is_bntl) && !tao->ls_jac) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualJacobianRoutine() must be called before setup!"); 230e1e80dc8SAlp Dener if (!tao->gradient) { 231e1e80dc8SAlp Dener ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 232e1e80dc8SAlp Dener } 233737f463aSAlp Dener if (!gn->x_work) { 234737f463aSAlp Dener ierr = VecDuplicate(tao->solution,&gn->x_work);CHKERRQ(ierr); 235737f463aSAlp Dener } 236737f463aSAlp Dener if (!gn->r_work) { 237737f463aSAlp Dener ierr = VecDuplicate(tao->ls_res,&gn->r_work);CHKERRQ(ierr); 238737f463aSAlp Dener } 239e1e80dc8SAlp Dener if (!gn->x_old) { 240e1e80dc8SAlp Dener ierr = VecDuplicate(tao->solution,&gn->x_old);CHKERRQ(ierr); 241e1e80dc8SAlp Dener ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr); 242e1e80dc8SAlp Dener } 2437cea06e1SXiang Huang 244470ec3f8SXiang Huang if (BRGN_REGULARIZATION_L1DICT == gn->reg_type) { 24530eeff36SXiang Huang if (gn->D) { 24630eeff36SXiang 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 */ 24730eeff36SXiang Huang } else { 24830eeff36SXiang Huang ierr = VecGetSize(tao->solution,&K);CHKERRQ(ierr); /* If user does not setup dict matrix, use identiy matrix, K=N */ 24930eeff36SXiang Huang } 2507cea06e1SXiang Huang if (!gn->y) { 2517cea06e1SXiang Huang ierr = VecCreate(PETSC_COMM_SELF,&gn->y);CHKERRQ(ierr); 2528e85b1b3SXiang Huang ierr = VecSetSizes(gn->y,PETSC_DECIDE,K);CHKERRQ(ierr); 2537cea06e1SXiang Huang ierr = VecSetFromOptions(gn->y);CHKERRQ(ierr); 2547cea06e1SXiang Huang ierr = VecSet(gn->y,0.0);CHKERRQ(ierr); 2557cea06e1SXiang Huang 2567cea06e1SXiang Huang } 2577cea06e1SXiang Huang if (!gn->y_work) { 2587cea06e1SXiang Huang ierr = VecDuplicate(gn->y,&gn->y_work);CHKERRQ(ierr); 2597cea06e1SXiang Huang } 2608ac80d48SXiang Huang if (!gn->diag) { 2617cea06e1SXiang Huang ierr = VecDuplicate(gn->y,&gn->diag);CHKERRQ(ierr); 2628ac80d48SXiang Huang ierr = VecSet(gn->diag,0.0);CHKERRQ(ierr); 2638ac80d48SXiang Huang } 26430eeff36SXiang Huang } 2650d71dc2bSXiang Huang 266e1e80dc8SAlp Dener if (!tao->setupcalled) { 267737f463aSAlp Dener /* Hessian setup */ 2688e85b1b3SXiang Huang ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 2698e85b1b3SXiang Huang ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 2708e85b1b3SXiang Huang ierr = MatSetSizes(gn->H,n,n,N,N);CHKERRQ(ierr); 271737f463aSAlp Dener ierr = MatSetType(gn->H,MATSHELL);CHKERRQ(ierr); 272737f463aSAlp Dener ierr = MatSetUp(gn->H);CHKERRQ(ierr); 273737f463aSAlp Dener ierr = MatShellSetOperation(gn->H,MATOP_MULT,(void (*)(void))GNHessianProd);CHKERRQ(ierr); 274737f463aSAlp Dener ierr = MatShellSetContext(gn->H,(void*)gn);CHKERRQ(ierr); 2758e85b1b3SXiang Huang /* Subsolver setup,include initial vector and dicttionary D */ 276e1e80dc8SAlp Dener ierr = TaoSetUpdate(gn->subsolver,GNHookFunction,(void*)gn);CHKERRQ(ierr); 277737f463aSAlp Dener ierr = TaoSetInitialVector(gn->subsolver,tao->solution);CHKERRQ(ierr); 278737f463aSAlp Dener if (tao->bounded) { 279737f463aSAlp Dener ierr = TaoSetVariableBounds(gn->subsolver,tao->XL,tao->XU);CHKERRQ(ierr); 280737f463aSAlp Dener } 281737f463aSAlp Dener ierr = TaoSetResidualRoutine(gn->subsolver,tao->ls_res,tao->ops->computeresidual,tao->user_lsresP);CHKERRQ(ierr); 2824ffbe8acSAlp Dener ierr = TaoSetJacobianResidualRoutine(gn->subsolver,tao->ls_jac,tao->ls_jac,tao->ops->computeresidualjacobian,tao->user_lsjacP);CHKERRQ(ierr); 283737f463aSAlp Dener ierr = TaoSetObjectiveAndGradientRoutine(gn->subsolver,GNObjectiveGradientEval,(void*)gn);CHKERRQ(ierr); 284737f463aSAlp Dener ierr = TaoSetHessianRoutine(gn->subsolver,gn->H,gn->H,GNComputeHessian,(void*)gn);CHKERRQ(ierr); 285e1e80dc8SAlp Dener /* Propagate some options down */ 286e1e80dc8SAlp Dener ierr = TaoSetTolerances(gn->subsolver,tao->gatol,tao->grtol,tao->gttol);CHKERRQ(ierr); 287e1e80dc8SAlp Dener ierr = TaoSetMaximumIterations(gn->subsolver,tao->max_it);CHKERRQ(ierr); 288e1e80dc8SAlp Dener ierr = TaoSetMaximumFunctionEvaluations(gn->subsolver,tao->max_funcs);CHKERRQ(ierr); 289737f463aSAlp Dener for (i=0; i<tao->numbermonitors; ++i) { 290737f463aSAlp Dener ierr = TaoSetMonitor(gn->subsolver,tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i]);CHKERRQ(ierr); 291737f463aSAlp Dener ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr); 292737f463aSAlp Dener } 293737f463aSAlp Dener ierr = TaoSetUp(gn->subsolver);CHKERRQ(ierr); 294e1e80dc8SAlp Dener } 295737f463aSAlp Dener PetscFunctionReturn(0); 296737f463aSAlp Dener } 297737f463aSAlp Dener 298737f463aSAlp Dener static PetscErrorCode TaoDestroy_BRGN(Tao tao) 299737f463aSAlp Dener { 300737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 301737f463aSAlp Dener PetscErrorCode ierr; 302737f463aSAlp Dener 303737f463aSAlp Dener PetscFunctionBegin; 304737f463aSAlp Dener if (tao->setupcalled) { 305e1e80dc8SAlp Dener ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr); 306737f463aSAlp Dener ierr = VecDestroy(&gn->x_work);CHKERRQ(ierr); 307737f463aSAlp Dener ierr = VecDestroy(&gn->r_work);CHKERRQ(ierr); 308e1e80dc8SAlp Dener ierr = VecDestroy(&gn->x_old);CHKERRQ(ierr); 3098ac80d48SXiang Huang ierr = VecDestroy(&gn->diag);CHKERRQ(ierr); 3107cea06e1SXiang Huang ierr = VecDestroy(&gn->y);CHKERRQ(ierr); 3117cea06e1SXiang Huang ierr = VecDestroy(&gn->y_work);CHKERRQ(ierr); 312737f463aSAlp Dener } 313737f463aSAlp Dener ierr = MatDestroy(&gn->H);CHKERRQ(ierr); 3147cea06e1SXiang Huang ierr = MatDestroy(&gn->D);CHKERRQ(ierr); 315463fc0ecSAlp Dener ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr); 316737f463aSAlp Dener ierr = TaoDestroy(&gn->subsolver);CHKERRQ(ierr); 317e1e80dc8SAlp Dener gn->parent = NULL; 318737f463aSAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 319737f463aSAlp Dener PetscFunctionReturn(0); 320737f463aSAlp Dener } 321737f463aSAlp Dener 3223850be85SAlp Dener /*MC 3233850be85SAlp Dener TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares 3243850be85SAlp Dener problems with bound constraints. This algorithm is a thin wrapper around TAOBNTL 325463fc0ecSAlp Dener that constructs the Gauss-Newton problem with the user-provided least-squares 326463fc0ecSAlp Dener residual and Jacobian. The algorithm offers both an L2-norm proximal point ("l2prox") 32701b716f5SXiang Huang regularizer, and a L1-norm dictionary regularizer ("l1dict"), where we approximate the 32801b716f5SXiang Huang L1-norm ||x||_1 by sum_i(sqrt(x_i^2+epsilon^2)-epsilon) with a small positive number epsilon. 32901b716f5SXiang Huang The user can also provide own regularization function. 3303850be85SAlp Dener 3313850be85SAlp Dener Options Database Keys: 33288fa4459SAlp Dener + -tao_brgn_regularizer_weight - regularizer weight (default 1e-4) 33388fa4459SAlp Dener . -tao_brgn_l1_smooth_epsilon - L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6) 33401b716f5SXiang Huang - -tao_brgn_regularization_type - regularization type ("user", "l2prox", "l1dict") 3353850be85SAlp Dener 3363850be85SAlp Dener Level: beginner 3373850be85SAlp Dener M*/ 338737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao) 339737f463aSAlp Dener { 340737f463aSAlp Dener TAO_BRGN *gn; 341737f463aSAlp Dener PetscErrorCode ierr; 342737f463aSAlp Dener 343737f463aSAlp Dener PetscFunctionBegin; 344737f463aSAlp Dener ierr = PetscNewLog(tao,&gn);CHKERRQ(ierr); 345737f463aSAlp Dener 346737f463aSAlp Dener tao->ops->destroy = TaoDestroy_BRGN; 347737f463aSAlp Dener tao->ops->setup = TaoSetUp_BRGN; 348737f463aSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BRGN; 349737f463aSAlp Dener tao->ops->view = TaoView_BRGN; 350737f463aSAlp Dener tao->ops->solve = TaoSolve_BRGN; 351737f463aSAlp Dener 352737f463aSAlp Dener tao->data = (void*)gn; 353e1e80dc8SAlp Dener gn->lambda = 1e-4; 3548ac80d48SXiang Huang gn->epsilon = 1e-6; 355e1e80dc8SAlp Dener gn->parent = tao; 356737f463aSAlp Dener 357737f463aSAlp Dener ierr = MatCreate(PetscObjectComm((PetscObject)tao),&gn->H);CHKERRQ(ierr); 358737f463aSAlp Dener ierr = MatSetOptionsPrefix(gn->H,"tao_brgn_hessian_");CHKERRQ(ierr); 359737f463aSAlp Dener 360737f463aSAlp Dener ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&gn->subsolver);CHKERRQ(ierr); 361737f463aSAlp Dener ierr = TaoSetType(gn->subsolver,TAOBNLS);CHKERRQ(ierr); 362737f463aSAlp Dener ierr = TaoSetOptionsPrefix(gn->subsolver,"tao_brgn_subsolver_");CHKERRQ(ierr); 363e1e80dc8SAlp Dener PetscFunctionReturn(0); 364e1e80dc8SAlp Dener } 365e1e80dc8SAlp Dener 366463fc0ecSAlp Dener /*@ 367e1e80dc8SAlp Dener TaoBRGNGetSubsolver - Get the pointer to the subsolver inside BRGN 368e1e80dc8SAlp Dener 369e1e80dc8SAlp Dener Collective on Tao 370e1e80dc8SAlp Dener 371463fc0ecSAlp Dener Level: advanced 372e1e80dc8SAlp Dener 373e1e80dc8SAlp Dener Input Parameters: 374e1e80dc8SAlp Dener + tao - the Tao solver context 375e1e80dc8SAlp Dener - subsolver - the Tao sub-solver context 376e1e80dc8SAlp Dener @*/ 377e1e80dc8SAlp Dener PetscErrorCode TaoBRGNGetSubsolver(Tao tao,Tao *subsolver) 378e1e80dc8SAlp Dener { 379e1e80dc8SAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 380e1e80dc8SAlp Dener 381e1e80dc8SAlp Dener PetscFunctionBegin; 382e1e80dc8SAlp Dener *subsolver = gn->subsolver; 383737f463aSAlp Dener PetscFunctionReturn(0); 384737f463aSAlp Dener } 385737f463aSAlp Dener 386463fc0ecSAlp Dener /*@ 387463fc0ecSAlp Dener TaoBRGNSetRegularizerWeight - Set the regularizer weight for the Gauss-Newton least-squares algorithm 388737f463aSAlp Dener 389737f463aSAlp Dener Collective on Tao 390737f463aSAlp Dener 391737f463aSAlp Dener Input Parameters: 392737f463aSAlp Dener + tao - the Tao solver context 3938e85b1b3SXiang Huang - lambda - L1-norm regularizer weight 394463fc0ecSAlp Dener 395463fc0ecSAlp Dener Level: beginner 396737f463aSAlp Dener @*/ 397a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao,PetscReal lambda) 398737f463aSAlp Dener { 399737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 400737f463aSAlp Dener 4018ac80d48SXiang Huang /* Initialize lambda here */ 4020d71dc2bSXiang Huang 403737f463aSAlp Dener PetscFunctionBegin; 404737f463aSAlp Dener gn->lambda = lambda; 405737f463aSAlp Dener PetscFunctionReturn(0); 406737f463aSAlp Dener } 4070d71dc2bSXiang Huang 408463fc0ecSAlp Dener /*@ 4098ac80d48SXiang Huang TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm 4108ac80d48SXiang Huang 4118ac80d48SXiang Huang Collective on Tao 4128ac80d48SXiang Huang 4138ac80d48SXiang Huang Input Parameters: 4148ac80d48SXiang Huang + tao - the Tao solver context 4158ac80d48SXiang Huang - epsilon - L1-norm smooth approximation parameter 416463fc0ecSAlp Dener 417463fc0ecSAlp Dener Level: advanced 4188ac80d48SXiang Huang @*/ 4198ac80d48SXiang Huang PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao,PetscReal epsilon) 4208ac80d48SXiang Huang { 4218ac80d48SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)tao->data; 4228ac80d48SXiang Huang 4238ac80d48SXiang Huang /* Initialize epsilon here */ 4248ac80d48SXiang Huang 4258ac80d48SXiang Huang PetscFunctionBegin; 4268ac80d48SXiang Huang gn->epsilon = epsilon; 4278ac80d48SXiang Huang PetscFunctionReturn(0); 4288ac80d48SXiang Huang } 4298e85b1b3SXiang Huang 430463fc0ecSAlp Dener /*@ 4318e85b1b3SXiang Huang TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem) 4328e85b1b3SXiang Huang 4338e85b1b3SXiang Huang Input Parameters: 4348e85b1b3SXiang Huang + tao - the Tao context 4351fc140a9SXiang Huang . dict - the user specified dictionary matrix. We allow to set a null dictionary, which means identity matrix by default 4368e85b1b3SXiang Huang 437463fc0ecSAlp Dener Level: advanced 4388e85b1b3SXiang Huang @*/ 4398e85b1b3SXiang Huang PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao,Mat dict) 4408e85b1b3SXiang Huang { 4418e85b1b3SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)tao->data; 4428e85b1b3SXiang Huang PetscErrorCode ierr; 4438e85b1b3SXiang Huang PetscFunctionBegin; 4448e85b1b3SXiang Huang PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 4458e85b1b3SXiang Huang if (dict) { 4468e85b1b3SXiang Huang PetscValidHeaderSpecific(dict,MAT_CLASSID,2); 447a3c390cfSAlp Dener PetscCheckSameComm(tao,1,dict,2); 4488e85b1b3SXiang Huang ierr = PetscObjectReference((PetscObject)dict);CHKERRQ(ierr); 4498e85b1b3SXiang Huang } 4508e85b1b3SXiang Huang ierr = MatDestroy(&gn->D);CHKERRQ(ierr); 4511fc140a9SXiang Huang gn->D = dict; 4528e85b1b3SXiang Huang PetscFunctionReturn(0); 4538e85b1b3SXiang Huang } 4548e85b1b3SXiang Huang 455a3c390cfSAlp Dener /*@C 456463fc0ecSAlp Dener TaoBRGNSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back 457463fc0ecSAlp Dener function into the algorithm. 458463fc0ecSAlp Dener 459463fc0ecSAlp Dener Input Parameters: 460463fc0ecSAlp Dener + tao - the Tao context 461463fc0ecSAlp Dener . func - function pointer for the regularizer value and gradient evaluation 462463fc0ecSAlp Dener - ctx - user context for the regularizer 463463fc0ecSAlp Dener 464463fc0ecSAlp Dener Level: advanced 465a3c390cfSAlp Dener @*/ 466a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (*func)(Tao,Vec,PetscReal *,Vec,void*),void *ctx) 467a3c390cfSAlp Dener { 468a3c390cfSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 469a3c390cfSAlp Dener 470a3c390cfSAlp Dener PetscFunctionBegin; 471a3c390cfSAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 472a3c390cfSAlp Dener if (ctx) { 473a3c390cfSAlp Dener gn->reg_obj_ctx = ctx; 474a3c390cfSAlp Dener } 475a3c390cfSAlp Dener if (func) { 476a3c390cfSAlp Dener gn->regularizerobjandgrad = func; 477a3c390cfSAlp Dener } 478a3c390cfSAlp Dener PetscFunctionReturn(0); 479a3c390cfSAlp Dener } 480a3c390cfSAlp Dener 481a3c390cfSAlp Dener /*@C 482463fc0ecSAlp Dener TaoBRGNSetRegularizerHessianRoutine - Sets the user-defined regularizer call-back 483463fc0ecSAlp Dener function into the algorithm. 484463fc0ecSAlp Dener 485463fc0ecSAlp Dener Input Parameters: 486463fc0ecSAlp Dener + tao - the Tao context 487463fc0ecSAlp Dener . Hreg - user-created matrix for the Hessian of the regularization term 488463fc0ecSAlp Dener . func - function pointer for the regularizer Hessian evaluation 489463fc0ecSAlp Dener - ctx - user context for the regularizer Hessian 490463fc0ecSAlp Dener 491463fc0ecSAlp Dener Level: advanced 492a3c390cfSAlp Dener @*/ 493a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao,Mat Hreg,PetscErrorCode (*func)(Tao,Vec,Mat,void*),void *ctx) 494a3c390cfSAlp Dener { 495a3c390cfSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 496a3c390cfSAlp Dener PetscErrorCode ierr; 497a3c390cfSAlp Dener 498a3c390cfSAlp Dener PetscFunctionBegin; 499a3c390cfSAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 500a3c390cfSAlp Dener if (Hreg) { 501a3c390cfSAlp Dener PetscValidHeaderSpecific(Hreg,MAT_CLASSID,2); 502a3c390cfSAlp Dener PetscCheckSameComm(tao,1,Hreg,2); 50301b716f5SXiang Huang } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONG,"NULL Hessian detected! User must provide valid Hessian for the regularizer."); 504a3c390cfSAlp Dener if (ctx) { 505a3c390cfSAlp Dener gn->reg_hess_ctx = ctx; 506a3c390cfSAlp Dener } 507a3c390cfSAlp Dener if (func) { 508a3c390cfSAlp Dener gn->regularizerhessian = func; 509a3c390cfSAlp Dener } 510a3c390cfSAlp Dener if (Hreg) { 511a3c390cfSAlp Dener ierr = PetscObjectReference((PetscObject)Hreg);CHKERRQ(ierr); 512a3c390cfSAlp Dener ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr); 513a3c390cfSAlp Dener gn->Hreg = Hreg; 514a3c390cfSAlp Dener } 515a3c390cfSAlp Dener PetscFunctionReturn(0); 516a3c390cfSAlp Dener } 517