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 7cd1c4666STristan Konolige #define BRGN_REGULARIZATION_LM 4 8cd1c4666STristan Konolige #define BRGN_REGULARIZATION_TYPES 5 9a3c390cfSAlp Dener 10cd1c4666STristan 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; 47cd1c4666STristan Konolige case BRGN_REGULARIZATION_LM: 48cd1c4666STristan Konolige ierr = VecPointwiseMult(gn->x_work,gn->damping,in);CHKERRQ(ierr); 49cd1c4666STristan Konolige ierr = VecAXPY(out,1,gn->x_work);CHKERRQ(ierr); 50cd1c4666STristan Konolige break; 51a3c390cfSAlp Dener } 520d71dc2bSXiang Huang PetscFunctionReturn(0); 530d71dc2bSXiang Huang } 54cd1c4666STristan Konolige static PetscErrorCode ComputeDamping(TAO_BRGN *gn) 55cd1c4666STristan Konolige { 56cd1c4666STristan Konolige const PetscScalar *diag_ary; 57cd1c4666STristan Konolige PetscScalar *damping_ary; 58cd1c4666STristan Konolige PetscInt i,n; 59cd1c4666STristan Konolige PetscErrorCode ierr; 60cd1c4666STristan Konolige 61cd1c4666STristan Konolige PetscFunctionBegin; 62cd1c4666STristan Konolige /* update damping */ 63cd1c4666STristan Konolige ierr = VecGetArray(gn->damping,&damping_ary);CHKERRQ(ierr); 64cd1c4666STristan Konolige ierr = VecGetArrayRead(gn->diag,&diag_ary);CHKERRQ(ierr); 65cd1c4666STristan Konolige ierr = VecGetLocalSize(gn->damping,&n);CHKERRQ(ierr); 66cd1c4666STristan Konolige for (i=0; i<n; i++) { 67cd1c4666STristan Konolige damping_ary[i] = PetscClipInterval(diag_ary[i],PETSC_SQRT_MACHINE_EPSILON,PetscSqrtReal(PETSC_MAX_REAL)); 68cd1c4666STristan Konolige } 69cd1c4666STristan Konolige ierr = VecScale(gn->damping,gn->lambda);CHKERRQ(ierr); 70cd1c4666STristan Konolige ierr = VecRestoreArray(gn->damping,&damping_ary);CHKERRQ(ierr); 71cd1c4666STristan Konolige ierr = VecRestoreArrayRead(gn->diag,&diag_ary);CHKERRQ(ierr); 72cd1c4666STristan Konolige PetscFunctionReturn(0); 73cd1c4666STristan Konolige } 74cd1c4666STristan Konolige 75cd1c4666STristan Konolige PetscErrorCode TaoBRGNGetDampingVector(Tao tao,Vec *d) 76cd1c4666STristan Konolige { 77cd1c4666STristan Konolige TAO_BRGN *gn = (TAO_BRGN *)tao->data; 78cd1c4666STristan Konolige 79cd1c4666STristan Konolige PetscFunctionBegin; 80cd1c4666STristan 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."); 81cd1c4666STristan Konolige *d = gn->damping; 82cd1c4666STristan Konolige PetscFunctionReturn(0); 83cd1c4666STristan 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; 153cd1c4666STristan Konolige PetscInt i,n,cstart,cend; 154cd1c4666STristan 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); 1595eb5f4d6SAlp Dener if (gn->mat_explicit) { 1605eb5f4d6SAlp Dener ierr = MatTransposeMatMult(tao->ls_jac, tao->ls_jac, MAT_REUSE_MATRIX, PETSC_DEFAULT, &gn->H);CHKERRQ(ierr); 1615eb5f4d6SAlp Dener } 1620d71dc2bSXiang Huang 163a3c390cfSAlp Dener switch (gn->reg_type) { 164470ec3f8SXiang Huang case BRGN_REGULARIZATION_USER: 165a3c390cfSAlp Dener ierr = (*gn->regularizerhessian)(tao,X,gn->Hreg,gn->reg_hess_ctx);CHKERRQ(ierr); 1665eb5f4d6SAlp Dener if (gn->mat_explicit) { 1675eb5f4d6SAlp Dener ierr = MatAXPY(gn->H, 1.0, gn->Hreg, DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1685eb5f4d6SAlp Dener } 169a3c390cfSAlp Dener break; 170a1c74439SHansol Suh case BRGN_REGULARIZATION_L2PURE: 1715eb5f4d6SAlp Dener if (gn->mat_explicit) { 1725eb5f4d6SAlp Dener ierr = MatShift(gn->H, gn->lambda);CHKERRQ(ierr); 1735eb5f4d6SAlp Dener } 174a1c74439SHansol Suh break; 175470ec3f8SXiang Huang case BRGN_REGULARIZATION_L2PROX: 1765eb5f4d6SAlp Dener if (gn->mat_explicit) { 1775eb5f4d6SAlp Dener ierr = MatShift(gn->H, gn->lambda);CHKERRQ(ierr); 1785eb5f4d6SAlp Dener } 179a3c390cfSAlp Dener break; 180470ec3f8SXiang Huang case BRGN_REGULARIZATION_L1DICT: 1817cea06e1SXiang 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 */ 182a3c390cfSAlp Dener if (gn->D) { 1837cea06e1SXiang Huang ierr = MatMult(gn->D,X,gn->y);CHKERRQ(ierr);/* y = D*x */ 184a3c390cfSAlp Dener } else { 185a3c390cfSAlp Dener ierr = VecCopy(X,gn->y);CHKERRQ(ierr); 186a3c390cfSAlp Dener } 1877cea06e1SXiang Huang ierr = VecPointwiseMult(gn->y_work,gn->y,gn->y);CHKERRQ(ierr); 1887cea06e1SXiang Huang ierr = VecShift(gn->y_work,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 1897cea06e1SXiang Huang ierr = VecCopy(gn->y_work,gn->diag);CHKERRQ(ierr); /* gn->diag = y.^2+epsilon^2 */ 1907cea06e1SXiang Huang ierr = VecSqrtAbs(gn->y_work);CHKERRQ(ierr); /* gn->y_work = sqrt(y.^2+epsilon^2) */ 1917cea06e1SXiang Huang ierr = VecPointwiseMult(gn->diag,gn->y_work,gn->diag);CHKERRQ(ierr);/* gn->diag = sqrt(y.^2+epsilon^2).^3 */ 1928ac80d48SXiang Huang ierr = VecReciprocal(gn->diag);CHKERRQ(ierr); 1938ac80d48SXiang Huang ierr = VecScale(gn->diag,gn->epsilon*gn->epsilon);CHKERRQ(ierr); 1945eb5f4d6SAlp Dener if (gn->mat_explicit) { 1955eb5f4d6SAlp Dener ierr = MatDiagonalSet(gn->H, gn->diag, ADD_VALUES);CHKERRQ(ierr); 1965eb5f4d6SAlp Dener } 197a3c390cfSAlp Dener break; 198cd1c4666STristan Konolige case BRGN_REGULARIZATION_LM: 199cd1c4666STristan Konolige /* compute diagonal of J^T J */ 200cd1c4666STristan Konolige ierr = MatGetSize(gn->parent->ls_jac,NULL,&n);CHKERRQ(ierr); 201cd1c4666STristan Konolige ierr = PetscMalloc1(n,&cnorms);CHKERRQ(ierr); 202cd1c4666STristan Konolige ierr = MatGetColumnNorms(gn->parent->ls_jac,NORM_2,cnorms);CHKERRQ(ierr); 203cd1c4666STristan Konolige ierr = MatGetOwnershipRangeColumn(gn->parent->ls_jac,&cstart,&cend);CHKERRQ(ierr); 204cd1c4666STristan Konolige ierr = VecGetArray(gn->diag,&diag_ary);CHKERRQ(ierr); 205cd1c4666STristan Konolige for (i = 0; i < cend-cstart; i++) { 206cd1c4666STristan Konolige diag_ary[i] = cnorms[cstart+i] * cnorms[cstart+i]; 207cd1c4666STristan Konolige } 208cd1c4666STristan Konolige ierr = VecRestoreArray(gn->diag,&diag_ary);CHKERRQ(ierr); 209cd1c4666STristan Konolige ierr = PetscFree(cnorms);CHKERRQ(ierr); 210cd1c4666STristan Konolige ierr = ComputeDamping(gn);CHKERRQ(ierr); 2115eb5f4d6SAlp Dener if (gn->mat_explicit) { 2125eb5f4d6SAlp Dener ierr = MatDiagonalSet(gn->H, gn->damping, ADD_VALUES);CHKERRQ(ierr); 2135eb5f4d6SAlp Dener } 214cd1c4666STristan Konolige break; 215a3c390cfSAlp Dener } 216e1e80dc8SAlp Dener PetscFunctionReturn(0); 217e1e80dc8SAlp Dener } 218e1e80dc8SAlp Dener 2198fcddce6SStefano Zampini static PetscErrorCode GNHookFunction(Tao tao,PetscInt iter, void *ctx) 220e1e80dc8SAlp Dener { 2218fcddce6SStefano Zampini TAO_BRGN *gn = (TAO_BRGN *)ctx; 222e1e80dc8SAlp Dener PetscErrorCode ierr; 223e1e80dc8SAlp Dener 224e1e80dc8SAlp Dener PetscFunctionBegin; 225e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */ 226e1e80dc8SAlp Dener gn->parent->nfuncs = tao->nfuncs; 227e1e80dc8SAlp Dener gn->parent->ngrads = tao->ngrads; 228e1e80dc8SAlp Dener gn->parent->nfuncgrads = tao->nfuncgrads; 229e1e80dc8SAlp Dener gn->parent->nhess = tao->nhess; 230e1e80dc8SAlp Dener gn->parent->niter = tao->niter; 231e1e80dc8SAlp Dener gn->parent->ksp_its = tao->ksp_its; 232e1e80dc8SAlp Dener gn->parent->ksp_tot_its = tao->ksp_tot_its; 233cd1c4666STristan Konolige gn->parent->fc = tao->fc; 234e1e80dc8SAlp Dener ierr = TaoGetConvergedReason(tao,&gn->parent->reason);CHKERRQ(ierr); 235e1e80dc8SAlp Dener /* Update the solution vectors */ 236e1e80dc8SAlp Dener if (iter == 0) { 237e1e80dc8SAlp Dener ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr); 238e1e80dc8SAlp Dener } else { 239e1e80dc8SAlp Dener ierr = VecCopy(tao->solution,gn->x_old);CHKERRQ(ierr); 240e1e80dc8SAlp Dener ierr = VecCopy(tao->solution,gn->parent->solution);CHKERRQ(ierr); 241e1e80dc8SAlp Dener } 242e1e80dc8SAlp Dener /* Update the gradient */ 243e1e80dc8SAlp Dener ierr = VecCopy(tao->gradient,gn->parent->gradient);CHKERRQ(ierr); 244cd1c4666STristan Konolige 245cd1c4666STristan Konolige /* Update damping parameter for LM */ 246cd1c4666STristan Konolige if (gn->reg_type == BRGN_REGULARIZATION_LM) { 247cd1c4666STristan Konolige if (iter > 0) { 248cd1c4666STristan Konolige if (gn->fc_old > tao->fc) { 249cd1c4666STristan Konolige gn->lambda = gn->lambda * gn->downhill_lambda_change; 250cd1c4666STristan Konolige } else { 251cd1c4666STristan Konolige /* uphill step */ 252cd1c4666STristan Konolige gn->lambda = gn->lambda * gn->uphill_lambda_change; 253cd1c4666STristan Konolige } 254cd1c4666STristan Konolige } 255cd1c4666STristan Konolige gn->fc_old = tao->fc; 256cd1c4666STristan Konolige } 257cd1c4666STristan Konolige 258e1e80dc8SAlp Dener /* Call general purpose update function */ 259e1e80dc8SAlp Dener if (gn->parent->ops->update) { 2608fcddce6SStefano Zampini ierr = (*gn->parent->ops->update)(gn->parent,gn->parent->niter,gn->parent->user_update);CHKERRQ(ierr); 261737f463aSAlp Dener } 262737f463aSAlp Dener PetscFunctionReturn(0); 263737f463aSAlp Dener } 264737f463aSAlp Dener 265737f463aSAlp Dener static PetscErrorCode TaoSolve_BRGN(Tao tao) 266737f463aSAlp Dener { 267737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 268737f463aSAlp Dener PetscErrorCode ierr; 269737f463aSAlp Dener 270737f463aSAlp Dener PetscFunctionBegin; 271737f463aSAlp Dener ierr = TaoSolve(gn->subsolver);CHKERRQ(ierr); 272e1e80dc8SAlp Dener /* Update basic tao information from the subsolver */ 273e1e80dc8SAlp Dener tao->nfuncs = gn->subsolver->nfuncs; 274e1e80dc8SAlp Dener tao->ngrads = gn->subsolver->ngrads; 275e1e80dc8SAlp Dener tao->nfuncgrads = gn->subsolver->nfuncgrads; 276e1e80dc8SAlp Dener tao->nhess = gn->subsolver->nhess; 277e1e80dc8SAlp Dener tao->niter = gn->subsolver->niter; 278e1e80dc8SAlp Dener tao->ksp_its = gn->subsolver->ksp_its; 279e1e80dc8SAlp Dener tao->ksp_tot_its = gn->subsolver->ksp_tot_its; 280e1e80dc8SAlp Dener ierr = TaoGetConvergedReason(gn->subsolver,&tao->reason);CHKERRQ(ierr); 281e1e80dc8SAlp Dener /* Update vectors */ 282e1e80dc8SAlp Dener ierr = VecCopy(gn->subsolver->solution,tao->solution);CHKERRQ(ierr); 283e1e80dc8SAlp Dener ierr = VecCopy(gn->subsolver->gradient,tao->gradient);CHKERRQ(ierr); 284737f463aSAlp Dener PetscFunctionReturn(0); 285737f463aSAlp Dener } 286737f463aSAlp Dener 287737f463aSAlp Dener static PetscErrorCode TaoSetFromOptions_BRGN(PetscOptionItems *PetscOptionsObject,Tao tao) 288737f463aSAlp Dener { 289737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 290cd1c4666STristan Konolige TaoLineSearch ls; 291737f463aSAlp Dener PetscErrorCode ierr; 292737f463aSAlp Dener 293737f463aSAlp Dener PetscFunctionBegin; 294a3c390cfSAlp 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); 2955eb5f4d6SAlp Dener ierr = PetscOptionsBool("-tao_brgn_mat_explicit","switches the Hessian construction to be an explicit matrix rather than MATSHELL","",gn->mat_explicit,&gn->mat_explicit,NULL);CHKERRQ(ierr); 29688fa4459SAlp Dener ierr = PetscOptionsReal("-tao_brgn_regularizer_weight","regularizer weight (default 1e-4)","",gn->lambda,&gn->lambda,NULL);CHKERRQ(ierr); 29788fa4459SAlp 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); 298*1e1ea65dSPierre Jolivet 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);CHKERRQ(ierr); 299*1e1ea65dSPierre Jolivet 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);CHKERRQ(ierr); 300470ec3f8SXiang 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); 301737f463aSAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 302cd1c4666STristan Konolige /* set unit line search direction as the default when using the lm regularizer */ 303cd1c4666STristan Konolige if (gn->reg_type == BRGN_REGULARIZATION_LM) { 304cd1c4666STristan Konolige ierr = TaoGetLineSearch(gn->subsolver,&ls);CHKERRQ(ierr); 305cd1c4666STristan Konolige ierr = TaoLineSearchSetType(ls,TAOLINESEARCHUNIT);CHKERRQ(ierr); 306cd1c4666STristan Konolige } 307737f463aSAlp Dener ierr = TaoSetFromOptions(gn->subsolver);CHKERRQ(ierr); 308737f463aSAlp Dener PetscFunctionReturn(0); 309737f463aSAlp Dener } 310737f463aSAlp Dener 311737f463aSAlp Dener static PetscErrorCode TaoView_BRGN(Tao tao,PetscViewer viewer) 312737f463aSAlp Dener { 313737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 314737f463aSAlp Dener PetscErrorCode ierr; 315737f463aSAlp Dener 316737f463aSAlp Dener PetscFunctionBegin; 317e1e80dc8SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 318737f463aSAlp Dener ierr = TaoView(gn->subsolver,viewer);CHKERRQ(ierr); 319e1e80dc8SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 320737f463aSAlp Dener PetscFunctionReturn(0); 321737f463aSAlp Dener } 322737f463aSAlp Dener 323737f463aSAlp Dener static PetscErrorCode TaoSetUp_BRGN(Tao tao) 324737f463aSAlp Dener { 325737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 326737f463aSAlp Dener PetscErrorCode ierr; 327737f463aSAlp Dener PetscBool is_bnls,is_bntr,is_bntl; 3288e85b1b3SXiang Huang PetscInt i,n,N,K; /* dict has size K*N*/ 329737f463aSAlp Dener 330737f463aSAlp Dener PetscFunctionBegin; 331737f463aSAlp Dener if (!tao->ls_res) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualRoutine() must be called before setup!"); 332737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNLS,&is_bnls);CHKERRQ(ierr); 333737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTR,&is_bntr);CHKERRQ(ierr); 334737f463aSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)gn->subsolver,TAOBNTL,&is_bntl);CHKERRQ(ierr); 335737f463aSAlp Dener if ((is_bnls || is_bntr || is_bntl) && !tao->ls_jac) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ORDER,"TaoSetResidualJacobianRoutine() must be called before setup!"); 336e1e80dc8SAlp Dener if (!tao->gradient) { 337e1e80dc8SAlp Dener ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 338e1e80dc8SAlp Dener } 339737f463aSAlp Dener if (!gn->x_work) { 340737f463aSAlp Dener ierr = VecDuplicate(tao->solution,&gn->x_work);CHKERRQ(ierr); 341737f463aSAlp Dener } 342737f463aSAlp Dener if (!gn->r_work) { 343737f463aSAlp Dener ierr = VecDuplicate(tao->ls_res,&gn->r_work);CHKERRQ(ierr); 344737f463aSAlp Dener } 345e1e80dc8SAlp Dener if (!gn->x_old) { 346e1e80dc8SAlp Dener ierr = VecDuplicate(tao->solution,&gn->x_old);CHKERRQ(ierr); 347e1e80dc8SAlp Dener ierr = VecSet(gn->x_old,0.0);CHKERRQ(ierr); 348e1e80dc8SAlp Dener } 3497cea06e1SXiang Huang 350470ec3f8SXiang Huang if (BRGN_REGULARIZATION_L1DICT == gn->reg_type) { 3512036730cSSajid Ali if (!gn->y) { 35230eeff36SXiang Huang if (gn->D) { 35330eeff36SXiang 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 */ 3542036730cSSajid Ali ierr = MatCreateVecs(gn->D,NULL,&gn->y);CHKERRQ(ierr); 35530eeff36SXiang Huang } else { 3562036730cSSajid Ali ierr = VecDuplicate(tao->solution,&gn->y);CHKERRQ(ierr); /* If user does not setup dict matrix, use identiy matrix, K=N */ 35730eeff36SXiang Huang } 3587cea06e1SXiang Huang ierr = VecSet(gn->y,0.0);CHKERRQ(ierr); 3597cea06e1SXiang Huang } 3607cea06e1SXiang Huang if (!gn->y_work) { 3617cea06e1SXiang Huang ierr = VecDuplicate(gn->y,&gn->y_work);CHKERRQ(ierr); 3627cea06e1SXiang Huang } 3638ac80d48SXiang Huang if (!gn->diag) { 3647cea06e1SXiang Huang ierr = VecDuplicate(gn->y,&gn->diag);CHKERRQ(ierr); 3658ac80d48SXiang Huang ierr = VecSet(gn->diag,0.0);CHKERRQ(ierr); 3668ac80d48SXiang Huang } 36730eeff36SXiang Huang } 368cd1c4666STristan Konolige if (BRGN_REGULARIZATION_LM == gn->reg_type) { 369cd1c4666STristan Konolige if (!gn->diag) { 3705eb5f4d6SAlp Dener ierr = MatCreateVecs(tao->ls_jac,&gn->diag,NULL);CHKERRQ(ierr); 371cd1c4666STristan Konolige } 372cd1c4666STristan Konolige if (!gn->damping) { 3735eb5f4d6SAlp Dener ierr = MatCreateVecs(tao->ls_jac,&gn->damping,NULL);CHKERRQ(ierr); 374cd1c4666STristan Konolige } 375cd1c4666STristan Konolige } 3760d71dc2bSXiang Huang 377e1e80dc8SAlp Dener if (!tao->setupcalled) { 378737f463aSAlp Dener /* Hessian setup */ 3795eb5f4d6SAlp Dener if (gn->mat_explicit) { 3805eb5f4d6SAlp Dener ierr = TaoComputeResidualJacobian(tao,tao->solution,tao->ls_jac,tao->ls_jac_pre);CHKERRQ(ierr); 381*1e1ea65dSPierre Jolivet ierr = MatTransposeMatMult(tao->ls_jac, tao->ls_jac, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &gn->H);CHKERRQ(ierr); 3825eb5f4d6SAlp Dener } else { 3838e85b1b3SXiang Huang ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 3848e85b1b3SXiang Huang ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 3855eb5f4d6SAlp Dener ierr = MatCreate(PetscObjectComm((PetscObject)tao),&gn->H);CHKERRQ(ierr); 3868e85b1b3SXiang Huang ierr = MatSetSizes(gn->H,n,n,N,N);CHKERRQ(ierr); 387737f463aSAlp Dener ierr = MatSetType(gn->H,MATSHELL);CHKERRQ(ierr); 388*1e1ea65dSPierre Jolivet ierr = MatSetOption(gn->H, MAT_SYMMETRIC, PETSC_TRUE);CHKERRQ(ierr); 389737f463aSAlp Dener ierr = MatShellSetOperation(gn->H,MATOP_MULT,(void (*)(void))GNHessianProd);CHKERRQ(ierr); 390737f463aSAlp Dener ierr = MatShellSetContext(gn->H,(void*)gn);CHKERRQ(ierr); 3915eb5f4d6SAlp Dener } 3925eb5f4d6SAlp Dener ierr = MatSetUp(gn->H);CHKERRQ(ierr); 3938e85b1b3SXiang Huang /* Subsolver setup,include initial vector and dicttionary D */ 394e1e80dc8SAlp Dener ierr = TaoSetUpdate(gn->subsolver,GNHookFunction,(void*)gn);CHKERRQ(ierr); 395737f463aSAlp Dener ierr = TaoSetInitialVector(gn->subsolver,tao->solution);CHKERRQ(ierr); 396737f463aSAlp Dener if (tao->bounded) { 397737f463aSAlp Dener ierr = TaoSetVariableBounds(gn->subsolver,tao->XL,tao->XU);CHKERRQ(ierr); 398737f463aSAlp Dener } 399737f463aSAlp Dener ierr = TaoSetResidualRoutine(gn->subsolver,tao->ls_res,tao->ops->computeresidual,tao->user_lsresP);CHKERRQ(ierr); 4004ffbe8acSAlp Dener ierr = TaoSetJacobianResidualRoutine(gn->subsolver,tao->ls_jac,tao->ls_jac,tao->ops->computeresidualjacobian,tao->user_lsjacP);CHKERRQ(ierr); 401737f463aSAlp Dener ierr = TaoSetObjectiveAndGradientRoutine(gn->subsolver,GNObjectiveGradientEval,(void*)gn);CHKERRQ(ierr); 402737f463aSAlp Dener ierr = TaoSetHessianRoutine(gn->subsolver,gn->H,gn->H,GNComputeHessian,(void*)gn);CHKERRQ(ierr); 403e1e80dc8SAlp Dener /* Propagate some options down */ 404e1e80dc8SAlp Dener ierr = TaoSetTolerances(gn->subsolver,tao->gatol,tao->grtol,tao->gttol);CHKERRQ(ierr); 405e1e80dc8SAlp Dener ierr = TaoSetMaximumIterations(gn->subsolver,tao->max_it);CHKERRQ(ierr); 406e1e80dc8SAlp Dener ierr = TaoSetMaximumFunctionEvaluations(gn->subsolver,tao->max_funcs);CHKERRQ(ierr); 407737f463aSAlp Dener for (i=0; i<tao->numbermonitors; ++i) { 408737f463aSAlp Dener ierr = TaoSetMonitor(gn->subsolver,tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i]);CHKERRQ(ierr); 409737f463aSAlp Dener ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr); 410737f463aSAlp Dener } 411737f463aSAlp Dener ierr = TaoSetUp(gn->subsolver);CHKERRQ(ierr); 412e1e80dc8SAlp Dener } 413737f463aSAlp Dener PetscFunctionReturn(0); 414737f463aSAlp Dener } 415737f463aSAlp Dener 416737f463aSAlp Dener static PetscErrorCode TaoDestroy_BRGN(Tao tao) 417737f463aSAlp Dener { 418737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 419737f463aSAlp Dener PetscErrorCode ierr; 420737f463aSAlp Dener 421737f463aSAlp Dener PetscFunctionBegin; 422737f463aSAlp Dener if (tao->setupcalled) { 423e1e80dc8SAlp Dener ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr); 424737f463aSAlp Dener ierr = VecDestroy(&gn->x_work);CHKERRQ(ierr); 425737f463aSAlp Dener ierr = VecDestroy(&gn->r_work);CHKERRQ(ierr); 426e1e80dc8SAlp Dener ierr = VecDestroy(&gn->x_old);CHKERRQ(ierr); 4278ac80d48SXiang Huang ierr = VecDestroy(&gn->diag);CHKERRQ(ierr); 4287cea06e1SXiang Huang ierr = VecDestroy(&gn->y);CHKERRQ(ierr); 4297cea06e1SXiang Huang ierr = VecDestroy(&gn->y_work);CHKERRQ(ierr); 430737f463aSAlp Dener } 431cd1c4666STristan Konolige ierr = VecDestroy(&gn->damping);CHKERRQ(ierr); 432cd1c4666STristan Konolige ierr = VecDestroy(&gn->diag);CHKERRQ(ierr); 433737f463aSAlp Dener ierr = MatDestroy(&gn->H);CHKERRQ(ierr); 4347cea06e1SXiang Huang ierr = MatDestroy(&gn->D);CHKERRQ(ierr); 435463fc0ecSAlp Dener ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr); 436737f463aSAlp Dener ierr = TaoDestroy(&gn->subsolver);CHKERRQ(ierr); 437e1e80dc8SAlp Dener gn->parent = NULL; 438737f463aSAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 439737f463aSAlp Dener PetscFunctionReturn(0); 440737f463aSAlp Dener } 441737f463aSAlp Dener 4423850be85SAlp Dener /*MC 4433850be85SAlp Dener TAOBRGN - Bounded Regularized Gauss-Newton method for solving nonlinear least-squares 4443850be85SAlp Dener problems with bound constraints. This algorithm is a thin wrapper around TAOBNTL 445463fc0ecSAlp Dener that constructs the Gauss-Newton problem with the user-provided least-squares 44660bb7533SHansol Suh residual and Jacobian. The algorithm offers an L2-norm ("l2pure"), L2-norm proximal point ("l2prox") 44760bb7533SHansol Suh regularizer, and L1-norm dictionary regularizer ("l1dict"), where we approximate the 44801b716f5SXiang Huang L1-norm ||x||_1 by sum_i(sqrt(x_i^2+epsilon^2)-epsilon) with a small positive number epsilon. 449cd1c4666STristan Konolige Also offered is the "lm" regularizer which uses a scaled diagonal of J^T J. 450cd1c4666STristan Konolige With the "lm" regularizer, BRGN is a Levenberg-Marquardt optimizer. 45101b716f5SXiang Huang The user can also provide own regularization function. 4523850be85SAlp Dener 4533850be85SAlp Dener Options Database Keys: 454cd1c4666STristan Konolige + -tao_brgn_regularization_type - regularization type ("user", "l2prox", "l2pure", "l1dict", "lm") (default "l2prox") 455c061e8e2SXiang Huang . -tao_brgn_regularizer_weight - regularizer weight (default 1e-4) 456c061e8e2SXiang Huang - -tao_brgn_l1_smooth_epsilon - L1-norm smooth approximation parameter: ||x||_1 = sum(sqrt(x.^2+epsilon^2)-epsilon) (default 1e-6) 4573850be85SAlp Dener 4583850be85SAlp Dener Level: beginner 4593850be85SAlp Dener M*/ 460737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BRGN(Tao tao) 461737f463aSAlp Dener { 462737f463aSAlp Dener TAO_BRGN *gn; 463737f463aSAlp Dener PetscErrorCode ierr; 464737f463aSAlp Dener 465737f463aSAlp Dener PetscFunctionBegin; 466737f463aSAlp Dener ierr = PetscNewLog(tao,&gn);CHKERRQ(ierr); 467737f463aSAlp Dener 468737f463aSAlp Dener tao->ops->destroy = TaoDestroy_BRGN; 469737f463aSAlp Dener tao->ops->setup = TaoSetUp_BRGN; 470737f463aSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BRGN; 471737f463aSAlp Dener tao->ops->view = TaoView_BRGN; 472737f463aSAlp Dener tao->ops->solve = TaoSolve_BRGN; 473737f463aSAlp Dener 474737f463aSAlp Dener tao->data = (void*)gn; 475d8bf7057SXiang Huang gn->reg_type = BRGN_REGULARIZATION_L2PROX; 476e1e80dc8SAlp Dener gn->lambda = 1e-4; 4778ac80d48SXiang Huang gn->epsilon = 1e-6; 478cd1c4666STristan Konolige gn->downhill_lambda_change = 1./5.; 479cd1c4666STristan Konolige gn->uphill_lambda_change = 1.5; 480e1e80dc8SAlp Dener gn->parent = tao; 481737f463aSAlp Dener 482737f463aSAlp Dener ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&gn->subsolver);CHKERRQ(ierr); 483737f463aSAlp Dener ierr = TaoSetType(gn->subsolver,TAOBNLS);CHKERRQ(ierr); 484737f463aSAlp Dener ierr = TaoSetOptionsPrefix(gn->subsolver,"tao_brgn_subsolver_");CHKERRQ(ierr); 485e1e80dc8SAlp Dener PetscFunctionReturn(0); 486e1e80dc8SAlp Dener } 487e1e80dc8SAlp Dener 488463fc0ecSAlp Dener /*@ 489e1e80dc8SAlp Dener TaoBRGNGetSubsolver - Get the pointer to the subsolver inside BRGN 490e1e80dc8SAlp Dener 491e1e80dc8SAlp Dener Collective on Tao 492e1e80dc8SAlp Dener 493463fc0ecSAlp Dener Level: advanced 494e1e80dc8SAlp Dener 495e1e80dc8SAlp Dener Input Parameters: 496e1e80dc8SAlp Dener + tao - the Tao solver context 497e1e80dc8SAlp Dener - subsolver - the Tao sub-solver context 498e1e80dc8SAlp Dener @*/ 499e1e80dc8SAlp Dener PetscErrorCode TaoBRGNGetSubsolver(Tao tao,Tao *subsolver) 500e1e80dc8SAlp Dener { 501e1e80dc8SAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 502e1e80dc8SAlp Dener 503e1e80dc8SAlp Dener PetscFunctionBegin; 504e1e80dc8SAlp Dener *subsolver = gn->subsolver; 505737f463aSAlp Dener PetscFunctionReturn(0); 506737f463aSAlp Dener } 507737f463aSAlp Dener 508463fc0ecSAlp Dener /*@ 509463fc0ecSAlp Dener TaoBRGNSetRegularizerWeight - Set the regularizer weight for the Gauss-Newton least-squares algorithm 510737f463aSAlp Dener 511737f463aSAlp Dener Collective on Tao 512737f463aSAlp Dener 513737f463aSAlp Dener Input Parameters: 514737f463aSAlp Dener + tao - the Tao solver context 5158e85b1b3SXiang Huang - lambda - L1-norm regularizer weight 516463fc0ecSAlp Dener 517463fc0ecSAlp Dener Level: beginner 518737f463aSAlp Dener @*/ 519a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerWeight(Tao tao,PetscReal lambda) 520737f463aSAlp Dener { 521737f463aSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 522737f463aSAlp Dener 5238ac80d48SXiang Huang /* Initialize lambda here */ 5240d71dc2bSXiang Huang 525737f463aSAlp Dener PetscFunctionBegin; 526737f463aSAlp Dener gn->lambda = lambda; 527737f463aSAlp Dener PetscFunctionReturn(0); 528737f463aSAlp Dener } 5290d71dc2bSXiang Huang 530463fc0ecSAlp Dener /*@ 5318ac80d48SXiang Huang TaoBRGNSetL1SmoothEpsilon - Set the L1-norm smooth approximation parameter for L1-regularized least-squares algorithm 5328ac80d48SXiang Huang 5338ac80d48SXiang Huang Collective on Tao 5348ac80d48SXiang Huang 5358ac80d48SXiang Huang Input Parameters: 5368ac80d48SXiang Huang + tao - the Tao solver context 5378ac80d48SXiang Huang - epsilon - L1-norm smooth approximation parameter 538463fc0ecSAlp Dener 539463fc0ecSAlp Dener Level: advanced 5408ac80d48SXiang Huang @*/ 5418ac80d48SXiang Huang PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao tao,PetscReal epsilon) 5428ac80d48SXiang Huang { 5438ac80d48SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)tao->data; 5448ac80d48SXiang Huang 5458ac80d48SXiang Huang /* Initialize epsilon here */ 5468ac80d48SXiang Huang 5478ac80d48SXiang Huang PetscFunctionBegin; 5488ac80d48SXiang Huang gn->epsilon = epsilon; 5498ac80d48SXiang Huang PetscFunctionReturn(0); 5508ac80d48SXiang Huang } 5518e85b1b3SXiang Huang 552463fc0ecSAlp Dener /*@ 5538e85b1b3SXiang Huang TaoBRGNSetDictionaryMatrix - bind the dictionary matrix from user application context to gn->D, for compressed sensing (with least-squares problem) 5548e85b1b3SXiang Huang 5558e85b1b3SXiang Huang Input Parameters: 5568e85b1b3SXiang Huang + tao - the Tao context 557a2b725a8SWilliam Gropp - dict - the user specified dictionary matrix. We allow to set a null dictionary, which means identity matrix by default 5588e85b1b3SXiang Huang 559463fc0ecSAlp Dener Level: advanced 5608e85b1b3SXiang Huang @*/ 5618e85b1b3SXiang Huang PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao tao,Mat dict) 5628e85b1b3SXiang Huang { 5638e85b1b3SXiang Huang TAO_BRGN *gn = (TAO_BRGN *)tao->data; 5648e85b1b3SXiang Huang PetscErrorCode ierr; 5658e85b1b3SXiang Huang PetscFunctionBegin; 5668e85b1b3SXiang Huang PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 5678e85b1b3SXiang Huang if (dict) { 5688e85b1b3SXiang Huang PetscValidHeaderSpecific(dict,MAT_CLASSID,2); 569a3c390cfSAlp Dener PetscCheckSameComm(tao,1,dict,2); 5708e85b1b3SXiang Huang ierr = PetscObjectReference((PetscObject)dict);CHKERRQ(ierr); 5718e85b1b3SXiang Huang } 5728e85b1b3SXiang Huang ierr = MatDestroy(&gn->D);CHKERRQ(ierr); 5731fc140a9SXiang Huang gn->D = dict; 5748e85b1b3SXiang Huang PetscFunctionReturn(0); 5758e85b1b3SXiang Huang } 5768e85b1b3SXiang Huang 577a3c390cfSAlp Dener /*@C 578463fc0ecSAlp Dener TaoBRGNSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back 579463fc0ecSAlp Dener function into the algorithm. 580463fc0ecSAlp Dener 581463fc0ecSAlp Dener Input Parameters: 582463fc0ecSAlp Dener + tao - the Tao context 583463fc0ecSAlp Dener . func - function pointer for the regularizer value and gradient evaluation 584463fc0ecSAlp Dener - ctx - user context for the regularizer 585463fc0ecSAlp Dener 586463fc0ecSAlp Dener Level: advanced 587a3c390cfSAlp Dener @*/ 588a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (*func)(Tao,Vec,PetscReal *,Vec,void*),void *ctx) 589a3c390cfSAlp Dener { 590a3c390cfSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 591a3c390cfSAlp Dener 592a3c390cfSAlp Dener PetscFunctionBegin; 593a3c390cfSAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 594a3c390cfSAlp Dener if (ctx) { 595a3c390cfSAlp Dener gn->reg_obj_ctx = ctx; 596a3c390cfSAlp Dener } 597a3c390cfSAlp Dener if (func) { 598a3c390cfSAlp Dener gn->regularizerobjandgrad = func; 599a3c390cfSAlp Dener } 600a3c390cfSAlp Dener PetscFunctionReturn(0); 601a3c390cfSAlp Dener } 602a3c390cfSAlp Dener 603a3c390cfSAlp Dener /*@C 604463fc0ecSAlp Dener TaoBRGNSetRegularizerHessianRoutine - Sets the user-defined regularizer call-back 605463fc0ecSAlp Dener function into the algorithm. 606463fc0ecSAlp Dener 607463fc0ecSAlp Dener Input Parameters: 608463fc0ecSAlp Dener + tao - the Tao context 609463fc0ecSAlp Dener . Hreg - user-created matrix for the Hessian of the regularization term 610463fc0ecSAlp Dener . func - function pointer for the regularizer Hessian evaluation 611463fc0ecSAlp Dener - ctx - user context for the regularizer Hessian 612463fc0ecSAlp Dener 613463fc0ecSAlp Dener Level: advanced 614a3c390cfSAlp Dener @*/ 615a3c390cfSAlp Dener PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao tao,Mat Hreg,PetscErrorCode (*func)(Tao,Vec,Mat,void*),void *ctx) 616a3c390cfSAlp Dener { 617a3c390cfSAlp Dener TAO_BRGN *gn = (TAO_BRGN *)tao->data; 618a3c390cfSAlp Dener PetscErrorCode ierr; 619a3c390cfSAlp Dener 620a3c390cfSAlp Dener PetscFunctionBegin; 621a3c390cfSAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 622a3c390cfSAlp Dener if (Hreg) { 623a3c390cfSAlp Dener PetscValidHeaderSpecific(Hreg,MAT_CLASSID,2); 624a3c390cfSAlp Dener PetscCheckSameComm(tao,1,Hreg,2); 62501b716f5SXiang Huang } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONG,"NULL Hessian detected! User must provide valid Hessian for the regularizer."); 626a3c390cfSAlp Dener if (ctx) { 627a3c390cfSAlp Dener gn->reg_hess_ctx = ctx; 628a3c390cfSAlp Dener } 629a3c390cfSAlp Dener if (func) { 630a3c390cfSAlp Dener gn->regularizerhessian = func; 631a3c390cfSAlp Dener } 632a3c390cfSAlp Dener if (Hreg) { 633a3c390cfSAlp Dener ierr = PetscObjectReference((PetscObject)Hreg);CHKERRQ(ierr); 634a3c390cfSAlp Dener ierr = MatDestroy(&gn->Hreg);CHKERRQ(ierr); 635a3c390cfSAlp Dener gn->Hreg = Hreg; 636a3c390cfSAlp Dener } 637a3c390cfSAlp Dener PetscFunctionReturn(0); 638a3c390cfSAlp Dener } 639