16285c0a3SHansol Suh #include <../src/tao/constrained/impls/admm/admm.h> /*I "petsctao.h" I*/ 26285c0a3SHansol Suh #include <petsctao.h> 36285c0a3SHansol Suh #include <petsc/private/petscimpl.h> 46285c0a3SHansol Suh 56285c0a3SHansol Suh /* Updates terminating criteria 66285c0a3SHansol Suh * 76285c0a3SHansol Suh * 1 ||r_k|| = ||Ax+Bz-c|| =< catol_admm* max{||Ax||,||Bz||,||c||} 86285c0a3SHansol Suh * 96285c0a3SHansol Suh * 2. Updates dual residual, d_k 106285c0a3SHansol Suh * 116285c0a3SHansol Suh * 3. ||d_k|| = ||mu*A^T*B(z_k-z_{k-1})|| =< gatol_admm * ||A^Ty|| */ 126285c0a3SHansol Suh 136285c0a3SHansol Suh static PetscBool cited = PETSC_FALSE; 146285c0a3SHansol Suh static const char citation[] = 156285c0a3SHansol Suh "@misc{xu2017adaptive,\n" 166285c0a3SHansol Suh " title={Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation},\n" 176285c0a3SHansol Suh " author={Zheng Xu and Mario A. T. Figueiredo and Xiaoming Yuan and Christoph Studer and Tom Goldstein},\n" 186285c0a3SHansol Suh " year={2017},\n" 196285c0a3SHansol Suh " eprint={1704.02712},\n" 206285c0a3SHansol Suh " archivePrefix={arXiv},\n" 216285c0a3SHansol Suh " primaryClass={cs.CV}\n" 226285c0a3SHansol Suh "} \n"; 236285c0a3SHansol Suh 24a82e8c82SStefano Zampini const char *const TaoADMMRegularizerTypes[] = {"REGULARIZER_USER","REGULARIZER_SOFT_THRESH","TaoADMMRegularizerType","TAO_ADMM_",NULL}; 25a82e8c82SStefano Zampini const char *const TaoADMMUpdateTypes[] = {"UPDATE_BASIC","UPDATE_ADAPTIVE","UPDATE_ADAPTIVE_RELAXED","TaoADMMUpdateType","TAO_ADMM_",NULL}; 26a82e8c82SStefano Zampini const char *const TaoALMMTypes[] = {"CLASSIC","PHR","TaoALMMType","TAO_ALMM_",NULL}; 27a82e8c82SStefano Zampini 286285c0a3SHansol Suh static PetscErrorCode TaoADMMToleranceUpdate(Tao tao) 296285c0a3SHansol Suh { 306285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 316285c0a3SHansol Suh PetscErrorCode ierr; 326285c0a3SHansol Suh PetscReal Axnorm,Bznorm,ATynorm,temp; 336285c0a3SHansol Suh Vec tempJR,tempL; 346285c0a3SHansol Suh Tao mis; 356285c0a3SHansol Suh 366285c0a3SHansol Suh PetscFunctionBegin; 376285c0a3SHansol Suh mis = am->subsolverX; 386285c0a3SHansol Suh tempJR = am->workJacobianRight; 396285c0a3SHansol Suh tempL = am->workLeft; 406285c0a3SHansol Suh /* ATy */ 416285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(mis, am->y, mis->jacobian_equality, mis->jacobian_equality_pre);CHKERRQ(ierr); 426285c0a3SHansol Suh ierr = MatMultTranspose(mis->jacobian_equality,am->y,tempJR);CHKERRQ(ierr); 436285c0a3SHansol Suh ierr = VecNorm(tempJR,NORM_2,&ATynorm);CHKERRQ(ierr); 446285c0a3SHansol Suh /* dualres = mu * ||AT(Bz-Bzold)||_2 */ 456285c0a3SHansol Suh ierr = VecWAXPY(tempJR,-1.,am->Bzold,am->Bz);CHKERRQ(ierr); 466285c0a3SHansol Suh ierr = MatMultTranspose(mis->jacobian_equality,tempJR,tempL);CHKERRQ(ierr); 476285c0a3SHansol Suh ierr = VecNorm(tempL,NORM_2,&am->dualres);CHKERRQ(ierr); 486285c0a3SHansol Suh am->dualres *= am->mu; 496285c0a3SHansol Suh 506285c0a3SHansol Suh /* ||Ax||_2, ||Bz||_2 */ 516285c0a3SHansol Suh ierr = VecNorm(am->Ax,NORM_2,&Axnorm);CHKERRQ(ierr); 526285c0a3SHansol Suh ierr = VecNorm(am->Bz,NORM_2,&Bznorm);CHKERRQ(ierr); 536285c0a3SHansol Suh 546285c0a3SHansol Suh /* Set catol to be catol_admm * max{||Ax||,||Bz||,||c||} * 556285c0a3SHansol Suh * Set gatol to be gatol_admm * ||A^Ty|| * 566285c0a3SHansol Suh * while cnorm is ||r_k||_2, and gnorm is ||d_k||_2 */ 576285c0a3SHansol Suh temp = am->catol_admm * PetscMax(Axnorm, (!am->const_norm) ? Bznorm : PetscMax(Bznorm,am->const_norm)); 586285c0a3SHansol Suh ierr = TaoSetConstraintTolerances(tao,temp,PETSC_DEFAULT);CHKERRQ(ierr); 596285c0a3SHansol Suh ierr = TaoSetTolerances(tao, am->gatol_admm*ATynorm, PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 606285c0a3SHansol Suh PetscFunctionReturn(0); 616285c0a3SHansol Suh } 626285c0a3SHansol Suh 636285c0a3SHansol Suh /* Penaly Update for Adaptive ADMM. */ 646285c0a3SHansol Suh static PetscErrorCode AdaptiveADMMPenaltyUpdate(Tao tao) 656285c0a3SHansol Suh { 666285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 676285c0a3SHansol Suh PetscErrorCode ierr; 686285c0a3SHansol Suh PetscReal ydiff_norm, yhatdiff_norm, Axdiff_norm, Bzdiff_norm, Axyhat, Bzy, a_sd, a_mg, a_k, b_sd, b_mg, b_k; 696285c0a3SHansol Suh PetscBool hflag, gflag; 706285c0a3SHansol Suh Vec tempJR,tempJR2; 716285c0a3SHansol Suh 726285c0a3SHansol Suh PetscFunctionBegin; 736285c0a3SHansol Suh tempJR = am->workJacobianRight; 746285c0a3SHansol Suh tempJR2 = am->workJacobianRight2; 756285c0a3SHansol Suh hflag = PETSC_FALSE; 766285c0a3SHansol Suh gflag = PETSC_FALSE; 776285c0a3SHansol Suh a_k = -1; 786285c0a3SHansol Suh b_k = -1; 796285c0a3SHansol Suh 806285c0a3SHansol Suh ierr = VecWAXPY(tempJR,-1.,am->Axold,am->Ax);CHKERRQ(ierr); 816285c0a3SHansol Suh ierr = VecWAXPY(tempJR2,-1.,am->yhatold,am->yhat);CHKERRQ(ierr); 826285c0a3SHansol Suh ierr = VecNorm(tempJR,NORM_2,&Axdiff_norm);CHKERRQ(ierr); 836285c0a3SHansol Suh ierr = VecNorm(tempJR2,NORM_2,&yhatdiff_norm);CHKERRQ(ierr); 846285c0a3SHansol Suh ierr = VecDot(tempJR,tempJR2,&Axyhat);CHKERRQ(ierr); 856285c0a3SHansol Suh 866285c0a3SHansol Suh ierr = VecWAXPY(tempJR,-1.,am->Bz0,am->Bz);CHKERRQ(ierr); 876285c0a3SHansol Suh ierr = VecWAXPY(tempJR2,-1.,am->y,am->y0);CHKERRQ(ierr); 886285c0a3SHansol Suh ierr = VecNorm(tempJR,NORM_2,&Bzdiff_norm);CHKERRQ(ierr); 896285c0a3SHansol Suh ierr = VecNorm(tempJR2,NORM_2,&ydiff_norm);CHKERRQ(ierr); 906285c0a3SHansol Suh ierr = VecDot(tempJR,tempJR2,&Bzy);CHKERRQ(ierr); 916285c0a3SHansol Suh 926285c0a3SHansol Suh if (Axyhat > am->orthval*Axdiff_norm*yhatdiff_norm + am->mueps) { 936285c0a3SHansol Suh hflag = PETSC_TRUE; 946285c0a3SHansol Suh a_sd = PetscSqr(yhatdiff_norm)/Axyhat; /* alphaSD */ 956285c0a3SHansol Suh a_mg = Axyhat/PetscSqr(Axdiff_norm); /* alphaMG */ 966285c0a3SHansol Suh a_k = (a_mg/a_sd) > 0.5 ? a_mg : a_sd - 0.5*a_mg; 976285c0a3SHansol Suh } 986285c0a3SHansol Suh if (Bzy > am->orthval*Bzdiff_norm*ydiff_norm + am->mueps) { 996285c0a3SHansol Suh gflag = PETSC_TRUE; 1006285c0a3SHansol Suh b_sd = PetscSqr(ydiff_norm)/Bzy; /* betaSD */ 1016285c0a3SHansol Suh b_mg = Bzy/PetscSqr(Bzdiff_norm); /* betaMG */ 1026285c0a3SHansol Suh b_k = (b_mg/b_sd) > 0.5 ? b_mg : b_sd - 0.5*b_mg; 1036285c0a3SHansol Suh } 1046285c0a3SHansol Suh am->muold = am->mu; 1056285c0a3SHansol Suh if (gflag && hflag) { 1066285c0a3SHansol Suh am->mu = PetscSqrtReal(a_k*b_k); 1076285c0a3SHansol Suh } else if (hflag) { 1086285c0a3SHansol Suh am->mu = a_k; 1096285c0a3SHansol Suh } else if (gflag) { 1106285c0a3SHansol Suh am->mu = b_k; 1116285c0a3SHansol Suh } 1126285c0a3SHansol Suh if (am->mu > am->muold) { 1136285c0a3SHansol Suh am->mu = am->muold; 1146285c0a3SHansol Suh } 1156285c0a3SHansol Suh if (am->mu < am->mumin) { 1166285c0a3SHansol Suh am->mu = am->mumin; 1176285c0a3SHansol Suh } 1186285c0a3SHansol Suh PetscFunctionReturn(0); 1196285c0a3SHansol Suh } 1206285c0a3SHansol Suh 1217f5c9be9SBarry Smith static PetscErrorCode TaoADMMSetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType type) 1226285c0a3SHansol Suh { 1236285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1246285c0a3SHansol Suh 1256285c0a3SHansol Suh PetscFunctionBegin; 1266285c0a3SHansol Suh am->regswitch = type; 1276285c0a3SHansol Suh PetscFunctionReturn(0); 1286285c0a3SHansol Suh } 1296285c0a3SHansol Suh 1307f5c9be9SBarry Smith static PetscErrorCode TaoADMMGetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType *type) 1316285c0a3SHansol Suh { 1326285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1336285c0a3SHansol Suh 1346285c0a3SHansol Suh PetscFunctionBegin; 1356285c0a3SHansol Suh *type = am->regswitch; 1366285c0a3SHansol Suh PetscFunctionReturn(0); 1376285c0a3SHansol Suh } 1386285c0a3SHansol Suh 1397f5c9be9SBarry Smith static PetscErrorCode TaoADMMSetUpdateType_ADMM(Tao tao, TaoADMMUpdateType type) 1406285c0a3SHansol Suh { 1416285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1426285c0a3SHansol Suh 1436285c0a3SHansol Suh PetscFunctionBegin; 1446285c0a3SHansol Suh am->update = type; 1456285c0a3SHansol Suh PetscFunctionReturn(0); 1466285c0a3SHansol Suh } 1476285c0a3SHansol Suh 1487f5c9be9SBarry Smith static PetscErrorCode TaoADMMGetUpdateType_ADMM(Tao tao, TaoADMMUpdateType *type) 1496285c0a3SHansol Suh { 1506285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1516285c0a3SHansol Suh 1526285c0a3SHansol Suh PetscFunctionBegin; 1536285c0a3SHansol Suh *type = am->update; 1546285c0a3SHansol Suh PetscFunctionReturn(0); 1556285c0a3SHansol Suh } 1566285c0a3SHansol Suh 1576285c0a3SHansol Suh /* This routine updates Jacobians with new x,z vectors, 1586285c0a3SHansol Suh * and then updates Ax and Bz vectors, then computes updated residual vector*/ 1596285c0a3SHansol Suh static PetscErrorCode ADMMUpdateConstraintResidualVector(Tao tao, Vec x, Vec z, Vec Ax, Vec Bz, Vec residual) 1606285c0a3SHansol Suh { 1616285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1626285c0a3SHansol Suh Tao mis,reg; 1636285c0a3SHansol Suh PetscErrorCode ierr; 1646285c0a3SHansol Suh 1656285c0a3SHansol Suh PetscFunctionBegin; 1666285c0a3SHansol Suh mis = am->subsolverX; 1676285c0a3SHansol Suh reg = am->subsolverZ; 1686285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(mis, x, mis->jacobian_equality, mis->jacobian_equality_pre);CHKERRQ(ierr); 1696285c0a3SHansol Suh ierr = MatMult(mis->jacobian_equality,x,Ax);CHKERRQ(ierr); 1706285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(reg, z, reg->jacobian_equality, reg->jacobian_equality_pre);CHKERRQ(ierr); 1716285c0a3SHansol Suh ierr = MatMult(reg->jacobian_equality,z,Bz);CHKERRQ(ierr); 1726285c0a3SHansol Suh 1736285c0a3SHansol Suh ierr = VecWAXPY(residual,1.,Bz,Ax);CHKERRQ(ierr); 1746285c0a3SHansol Suh if (am->constraint != NULL) { 1756285c0a3SHansol Suh ierr = VecAXPY(residual,-1.,am->constraint);CHKERRQ(ierr); 1766285c0a3SHansol Suh } 1776285c0a3SHansol Suh PetscFunctionReturn(0); 1786285c0a3SHansol Suh } 1796285c0a3SHansol Suh 1806285c0a3SHansol Suh /* Updates Augmented Lagrangians to given routines * 1816285c0a3SHansol Suh * For subsolverX, routine needs to be ComputeObjectiveAndGraidnet 1826285c0a3SHansol Suh * Separate Objective and Gradient routines are not supported. */ 1836285c0a3SHansol Suh static PetscErrorCode SubObjGradUpdate(Tao tao, Vec x, PetscReal *f, Vec g, void *ptr) 1846285c0a3SHansol Suh { 1856285c0a3SHansol Suh Tao parent = (Tao)ptr; 1866285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 1876285c0a3SHansol Suh PetscErrorCode ierr; 1886285c0a3SHansol Suh PetscReal temp,temp2; 1896285c0a3SHansol Suh Vec tempJR; 1906285c0a3SHansol Suh 1916285c0a3SHansol Suh PetscFunctionBegin; 1926285c0a3SHansol Suh tempJR = am->workJacobianRight; 1936285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(parent, x, am->subsolverZ->solution, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 1946285c0a3SHansol Suh ierr = (*am->ops->misfitobjgrad)(am->subsolverX,x,f,g,am->misfitobjgradP);CHKERRQ(ierr); 1956285c0a3SHansol Suh 1966285c0a3SHansol Suh am->last_misfit_val = *f; 1976285c0a3SHansol Suh /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 1986285c0a3SHansol Suh ierr = VecTDot(am->residual,am->y,&temp);CHKERRQ(ierr); 1996285c0a3SHansol Suh ierr = VecTDot(am->residual,am->residual,&temp2);CHKERRQ(ierr); 2006285c0a3SHansol Suh *f += temp + (am->mu/2)*temp2; 2016285c0a3SHansol Suh 2026285c0a3SHansol Suh /* Gradient. Add + mu*AT(Ax+Bz-c) + yTA*/ 2036285c0a3SHansol Suh ierr = MatMultTranspose(tao->jacobian_equality,am->residual,tempJR);CHKERRQ(ierr); 2046285c0a3SHansol Suh ierr = VecAXPY(g,am->mu,tempJR);CHKERRQ(ierr); 2056285c0a3SHansol Suh ierr = MatMultTranspose(tao->jacobian_equality,am->y,tempJR);CHKERRQ(ierr); 2066285c0a3SHansol Suh ierr = VecAXPY(g,1.,tempJR);CHKERRQ(ierr); 2076285c0a3SHansol Suh PetscFunctionReturn(0); 2086285c0a3SHansol Suh } 2096285c0a3SHansol Suh 2106285c0a3SHansol Suh /* Updates Augmented Lagrangians to given routines 2116285c0a3SHansol Suh * For subsolverZ, routine needs to be ComputeObjectiveAndGraidnet 2126285c0a3SHansol Suh * Separate Objective and Gradient routines are not supported. */ 2136285c0a3SHansol Suh static PetscErrorCode RegObjGradUpdate(Tao tao, Vec z, PetscReal *f, Vec g, void *ptr) 2146285c0a3SHansol Suh { 2156285c0a3SHansol Suh Tao parent = (Tao)ptr; 2166285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 2176285c0a3SHansol Suh PetscErrorCode ierr; 2186285c0a3SHansol Suh PetscReal temp,temp2; 2196285c0a3SHansol Suh Vec tempJR; 2206285c0a3SHansol Suh 2216285c0a3SHansol Suh PetscFunctionBegin; 2226285c0a3SHansol Suh tempJR = am->workJacobianRight; 2236285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(parent, am->subsolverX->solution, z, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 2246285c0a3SHansol Suh ierr = (*am->ops->regobjgrad)(am->subsolverZ,z,f,g,am->regobjgradP);CHKERRQ(ierr); 2256285c0a3SHansol Suh am->last_reg_val= *f; 2266285c0a3SHansol Suh /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 2276285c0a3SHansol Suh ierr = VecTDot(am->residual,am->y,&temp);CHKERRQ(ierr); 2286285c0a3SHansol Suh ierr = VecTDot(am->residual,am->residual,&temp2);CHKERRQ(ierr); 2296285c0a3SHansol Suh *f += temp + (am->mu/2)*temp2; 2306285c0a3SHansol Suh 2316285c0a3SHansol Suh /* Gradient. Add + mu*BT(Ax+Bz-c) + yTB*/ 2326285c0a3SHansol Suh ierr = MatMultTranspose(am->subsolverZ->jacobian_equality,am->residual,tempJR);CHKERRQ(ierr); 2336285c0a3SHansol Suh ierr = VecAXPY(g,am->mu,tempJR);CHKERRQ(ierr); 2346285c0a3SHansol Suh ierr = MatMultTranspose(am->subsolverZ->jacobian_equality,am->y,tempJR);CHKERRQ(ierr); 2356285c0a3SHansol Suh ierr = VecAXPY(g,1.,tempJR);CHKERRQ(ierr); 2366285c0a3SHansol Suh PetscFunctionReturn(0); 2376285c0a3SHansol Suh } 2386285c0a3SHansol Suh 2396285c0a3SHansol Suh /* Computes epsilon padded L1 norm lambda*sum(sqrt(x^2+eps^2)-eps */ 2406285c0a3SHansol Suh static PetscErrorCode ADMML1EpsilonNorm(Tao tao, Vec x, PetscReal eps, PetscReal *norm) 2416285c0a3SHansol Suh { 2426285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 2436285c0a3SHansol Suh PetscInt N; 2446285c0a3SHansol Suh PetscErrorCode ierr; 2456285c0a3SHansol Suh 2466285c0a3SHansol Suh PetscFunctionBegin; 2476285c0a3SHansol Suh ierr = VecGetSize(am->workLeft,&N);CHKERRQ(ierr); 2486285c0a3SHansol Suh ierr = VecPointwiseMult(am->workLeft,x,x);CHKERRQ(ierr); 2496285c0a3SHansol Suh ierr = VecShift(am->workLeft,am->l1epsilon*am->l1epsilon);CHKERRQ(ierr); 2506285c0a3SHansol Suh ierr = VecSqrtAbs(am->workLeft);CHKERRQ(ierr); 2516285c0a3SHansol Suh ierr = VecSum(am->workLeft,norm);CHKERRQ(ierr); 2526285c0a3SHansol Suh *norm += N*am->l1epsilon; 2536285c0a3SHansol Suh *norm *= am->lambda; 2546285c0a3SHansol Suh PetscFunctionReturn(0); 2556285c0a3SHansol Suh } 2566285c0a3SHansol Suh 2576285c0a3SHansol Suh static PetscErrorCode ADMMInternalHessianUpdate(Mat H, Mat Constraint, PetscBool Identity, void *ptr) 2586285c0a3SHansol Suh { 2596285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)ptr; 2606285c0a3SHansol Suh PetscErrorCode ierr; 2616285c0a3SHansol Suh 2626285c0a3SHansol Suh PetscFunctionBegin; 2636285c0a3SHansol Suh switch (am->update) { 2646285c0a3SHansol Suh case (TAO_ADMM_UPDATE_BASIC): 2656285c0a3SHansol Suh break; 2666285c0a3SHansol Suh case (TAO_ADMM_UPDATE_ADAPTIVE): 2676285c0a3SHansol Suh case (TAO_ADMM_UPDATE_ADAPTIVE_RELAXED): 2686285c0a3SHansol Suh if (H && (am->muold != am->mu)) { 2696285c0a3SHansol Suh if (!Identity) { 2706285c0a3SHansol Suh ierr = MatAXPY(H,am->mu-am->muold,Constraint,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2716285c0a3SHansol Suh } else { 2726285c0a3SHansol Suh ierr = MatShift(H,am->mu-am->muold);CHKERRQ(ierr); 2736285c0a3SHansol Suh } 2746285c0a3SHansol Suh } 2756285c0a3SHansol Suh break; 2766285c0a3SHansol Suh } 2776285c0a3SHansol Suh PetscFunctionReturn(0); 2786285c0a3SHansol Suh } 2796285c0a3SHansol Suh 2806285c0a3SHansol Suh /* Updates Hessian - adds second derivative of augmented Lagrangian 2816285c0a3SHansol Suh * H \gets H + \rho*ATA 2826285c0a3SHansol Suh * Here, \rho does not change in TAO_ADMM_UPDATE_BASIC - thus no-op 2836285c0a3SHansol Suh * For ADAPTAIVE,ADAPTIVE_RELAXED, 2846285c0a3SHansol Suh * H \gets H + (\rho-\rhoold)*ATA 2856285c0a3SHansol Suh * Here, we assume that A is linear constraint i.e., doesnt change. 2866285c0a3SHansol Suh * Thus, for both ADAPTIVE, and RELAXED, ATA matrix is pre-set (except for A=I (null case)) see TaoSetUp_ADMM */ 2876285c0a3SHansol Suh static PetscErrorCode SubHessianUpdate(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 2886285c0a3SHansol Suh { 2896285c0a3SHansol Suh Tao parent = (Tao)ptr; 2906285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 2916285c0a3SHansol Suh PetscErrorCode ierr; 2926285c0a3SHansol Suh 2936285c0a3SHansol Suh PetscFunctionBegin; 2946285c0a3SHansol Suh if (am->Hxchange) { 2956285c0a3SHansol Suh /* Case where Hessian gets updated with respect to x vector input. */ 2966285c0a3SHansol Suh ierr = (*am->ops->misfithess)(am->subsolverX,x,H,Hpre,am->misfithessP);CHKERRQ(ierr); 2976285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am);CHKERRQ(ierr); 2986285c0a3SHansol Suh } else if (am->Hxbool) { 2996285c0a3SHansol Suh /* Hessian doesn't get updated. H(x) = c */ 3006285c0a3SHansol Suh /* Update Lagrangian only only per TAO call */ 3016285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am);CHKERRQ(ierr); 3026285c0a3SHansol Suh am->Hxbool = PETSC_FALSE; 3036285c0a3SHansol Suh } 3046285c0a3SHansol Suh PetscFunctionReturn(0); 3056285c0a3SHansol Suh } 3066285c0a3SHansol Suh 3076285c0a3SHansol Suh /* Same as SubHessianUpdate, except for B matrix instead of A matrix */ 3086285c0a3SHansol Suh static PetscErrorCode RegHessianUpdate(Tao tao, Vec z, Mat H, Mat Hpre, void *ptr) 3096285c0a3SHansol Suh { 3106285c0a3SHansol Suh Tao parent = (Tao)ptr; 3116285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 3126285c0a3SHansol Suh PetscErrorCode ierr; 3136285c0a3SHansol Suh 3146285c0a3SHansol Suh PetscFunctionBegin; 3156285c0a3SHansol Suh 3166285c0a3SHansol Suh if (am->Hzchange) { 3176285c0a3SHansol Suh /* Case where Hessian gets updated with respect to x vector input. */ 3186285c0a3SHansol Suh ierr = (*am->ops->reghess)(am->subsolverZ,z,H,Hpre,am->reghessP);CHKERRQ(ierr); 3196285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am);CHKERRQ(ierr); 3206285c0a3SHansol Suh } else if (am->Hzbool) { 3216285c0a3SHansol Suh /* Hessian doesn't get updated. H(x) = c */ 3226285c0a3SHansol Suh /* Update Lagrangian only only per TAO call */ 3236285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am);CHKERRQ(ierr); 3246285c0a3SHansol Suh am->Hzbool = PETSC_FALSE; 3256285c0a3SHansol Suh } 3266285c0a3SHansol Suh PetscFunctionReturn(0); 3276285c0a3SHansol Suh } 3286285c0a3SHansol Suh 3296285c0a3SHansol Suh /* Shell Matrix routine for A matrix. 3306285c0a3SHansol Suh * This gets used when user puts NULL for 3316285c0a3SHansol Suh * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...) 3326285c0a3SHansol Suh * Essentially sets A=I*/ 3336285c0a3SHansol Suh static PetscErrorCode JacobianIdentity(Mat mat,Vec in,Vec out) 3346285c0a3SHansol Suh { 3356285c0a3SHansol Suh PetscErrorCode ierr; 3366285c0a3SHansol Suh 3376285c0a3SHansol Suh PetscFunctionBegin; 3386285c0a3SHansol Suh ierr = VecCopy(in,out);CHKERRQ(ierr); 3396285c0a3SHansol Suh PetscFunctionReturn(0); 3406285c0a3SHansol Suh } 3416285c0a3SHansol Suh 3426285c0a3SHansol Suh /* Shell Matrix routine for B matrix. 3436285c0a3SHansol Suh * This gets used when user puts NULL for 3446285c0a3SHansol Suh * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...) 3456285c0a3SHansol Suh * Sets B=-I */ 3466285c0a3SHansol Suh static PetscErrorCode JacobianIdentityB(Mat mat,Vec in,Vec out) 3476285c0a3SHansol Suh { 3486285c0a3SHansol Suh PetscErrorCode ierr; 3496285c0a3SHansol Suh 3506285c0a3SHansol Suh PetscFunctionBegin; 3516285c0a3SHansol Suh ierr = VecCopy(in,out);CHKERRQ(ierr); 3526285c0a3SHansol Suh ierr = VecScale(out,-1.);CHKERRQ(ierr); 3536285c0a3SHansol Suh PetscFunctionReturn(0); 3546285c0a3SHansol Suh } 3556285c0a3SHansol Suh 3566285c0a3SHansol Suh /* Solve f(x) + g(z) s.t. Ax + Bz = c */ 3576285c0a3SHansol Suh static PetscErrorCode TaoSolve_ADMM(Tao tao) 3586285c0a3SHansol Suh { 3596285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 3606285c0a3SHansol Suh PetscErrorCode ierr; 3616285c0a3SHansol Suh PetscInt N; 3626285c0a3SHansol Suh PetscReal reg_func; 3636285c0a3SHansol Suh PetscBool is_reg_shell; 3646285c0a3SHansol Suh Vec tempL; 3656285c0a3SHansol Suh 3666285c0a3SHansol Suh PetscFunctionBegin; 3676285c0a3SHansol Suh if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 368*3c859ba3SBarry Smith PetscCheck(am->subsolverX->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetMisfitConstraintJacobian() first"); 369*3c859ba3SBarry Smith PetscCheck(am->subsolverZ->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerConstraintJacobian() first"); 3706285c0a3SHansol Suh if (am->constraint != NULL) { 3716285c0a3SHansol Suh ierr = VecNorm(am->constraint,NORM_2,&am->const_norm);CHKERRQ(ierr); 3726285c0a3SHansol Suh } 3736285c0a3SHansol Suh } 3746285c0a3SHansol Suh tempL = am->workLeft; 3756285c0a3SHansol Suh ierr = VecGetSize(tempL,&N);CHKERRQ(ierr); 3766285c0a3SHansol Suh 3776285c0a3SHansol Suh if (am->Hx && am->ops->misfithess) { 378a82e8c82SStefano Zampini ierr = TaoSetHessian(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao);CHKERRQ(ierr); 3796285c0a3SHansol Suh } 3806285c0a3SHansol Suh 3816285c0a3SHansol Suh if (!am->zJI) { 3826285c0a3SHansol Suh /* Currently, B is assumed to be a linear system, i.e., not getting updated*/ 3836285c0a3SHansol Suh ierr = MatTransposeMatMult(am->JB,am->JB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->BTB));CHKERRQ(ierr); 3846285c0a3SHansol Suh } 3856285c0a3SHansol Suh if (!am->xJI) { 3866285c0a3SHansol Suh /* Currently, A is assumed to be a linear system, i.e., not getting updated*/ 3876285c0a3SHansol Suh ierr = MatTransposeMatMult(am->subsolverX->jacobian_equality,am->subsolverX->jacobian_equality,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->ATA));CHKERRQ(ierr); 3886285c0a3SHansol Suh } 3896285c0a3SHansol Suh 3906285c0a3SHansol Suh is_reg_shell = PETSC_FALSE; 3916285c0a3SHansol Suh 3926285c0a3SHansol Suh ierr = PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell);CHKERRQ(ierr); 3936285c0a3SHansol Suh 3946285c0a3SHansol Suh if (!is_reg_shell) { 3956285c0a3SHansol Suh switch (am->regswitch) { 3966285c0a3SHansol Suh case (TAO_ADMM_REGULARIZER_USER): 397*3c859ba3SBarry Smith PetscCheck(am->ops->regobjgrad,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerObjectiveAndGradientRoutine() first if one wishes to use TAO_ADMM_REGULARIZER_USER with non-TAOSHELL type"); 3986285c0a3SHansol Suh break; 3996285c0a3SHansol Suh case (TAO_ADMM_REGULARIZER_SOFT_THRESH): 4006285c0a3SHansol Suh /* Soft Threshold. */ 4016285c0a3SHansol Suh break; 4026285c0a3SHansol Suh } 4036285c0a3SHansol Suh if (am->ops->regobjgrad) { 404a82e8c82SStefano Zampini ierr = TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao);CHKERRQ(ierr); 4056285c0a3SHansol Suh } 4066285c0a3SHansol Suh if (am->Hz && am->ops->reghess) { 407a82e8c82SStefano Zampini ierr = TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao);CHKERRQ(ierr); 4086285c0a3SHansol Suh } 4096285c0a3SHansol Suh } 4106285c0a3SHansol Suh 4116285c0a3SHansol Suh switch (am->update) { 4126285c0a3SHansol Suh case TAO_ADMM_UPDATE_BASIC: 4136285c0a3SHansol Suh if (am->subsolverX->hessian) { 4146285c0a3SHansol Suh /* In basic case, Hessian does not get updated w.r.t. to spectral penalty 4156285c0a3SHansol Suh * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/ 4166285c0a3SHansol Suh if (!am->xJI) { 4176285c0a3SHansol Suh ierr = MatAXPY(am->subsolverX->hessian,am->mu,am->ATA,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 4186285c0a3SHansol Suh } else { 4196285c0a3SHansol Suh ierr = MatShift(am->subsolverX->hessian,am->mu);CHKERRQ(ierr); 4206285c0a3SHansol Suh } 4216285c0a3SHansol Suh } 4226285c0a3SHansol Suh if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) { 4236285c0a3SHansol Suh if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) { 4246285c0a3SHansol Suh ierr = MatAXPY(am->subsolverZ->hessian,am->mu,am->BTB,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 4256285c0a3SHansol Suh } else { 4266285c0a3SHansol Suh ierr = MatShift(am->subsolverZ->hessian,am->mu);CHKERRQ(ierr); 4276285c0a3SHansol Suh } 4286285c0a3SHansol Suh } 4296285c0a3SHansol Suh break; 4306285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE: 4316285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 4326285c0a3SHansol Suh break; 4336285c0a3SHansol Suh } 4346285c0a3SHansol Suh 4356285c0a3SHansol Suh ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 4366285c0a3SHansol Suh tao->reason = TAO_CONTINUE_ITERATING; 4376285c0a3SHansol Suh 4386285c0a3SHansol Suh while (tao->reason == TAO_CONTINUE_ITERATING) { 4396285c0a3SHansol Suh if (tao->ops->update) { 4406285c0a3SHansol Suh ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 4416285c0a3SHansol Suh } 4426285c0a3SHansol Suh ierr = VecCopy(am->Bz, am->Bzold);CHKERRQ(ierr); 4436285c0a3SHansol Suh 4446285c0a3SHansol Suh /* x update */ 4456285c0a3SHansol Suh ierr = TaoSolve(am->subsolverX);CHKERRQ(ierr); 4466285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre);CHKERRQ(ierr); 4476285c0a3SHansol Suh ierr = MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution,am->Ax);CHKERRQ(ierr); 4486285c0a3SHansol Suh 4496285c0a3SHansol Suh am->Hxbool = PETSC_TRUE; 4506285c0a3SHansol Suh 4516285c0a3SHansol Suh /* z update */ 4526285c0a3SHansol Suh switch (am->regswitch) { 4536285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_USER: 4546285c0a3SHansol Suh ierr = TaoSolve(am->subsolverZ);CHKERRQ(ierr); 4556285c0a3SHansol Suh break; 4566285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_SOFT_THRESH: 4576285c0a3SHansol Suh /* L1 assumes A,B jacobians are identity nxn matrix */ 4586285c0a3SHansol Suh ierr = VecWAXPY(am->workJacobianRight,1/am->mu,am->y,am->Ax);CHKERRQ(ierr); 4596285c0a3SHansol Suh ierr = TaoSoftThreshold(am->workJacobianRight,-am->lambda/am->mu,am->lambda/am->mu,am->subsolverZ->solution);CHKERRQ(ierr); 4606285c0a3SHansol Suh break; 4616285c0a3SHansol Suh } 4626285c0a3SHansol Suh am->Hzbool = PETSC_TRUE; 4636285c0a3SHansol Suh /* Returns Ax + Bz - c with updated Ax,Bz vectors */ 4646285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 4656285c0a3SHansol Suh /* Dual variable, y += y + mu*(Ax+Bz-c) */ 4666285c0a3SHansol Suh ierr = VecWAXPY(am->y, am->mu, am->residual, am->yold);CHKERRQ(ierr); 4676285c0a3SHansol Suh 4686285c0a3SHansol Suh /* stopping tolerance update */ 4696285c0a3SHansol Suh ierr = TaoADMMToleranceUpdate(tao);CHKERRQ(ierr); 4706285c0a3SHansol Suh 4716285c0a3SHansol Suh /* Updating Spectral Penalty */ 4726285c0a3SHansol Suh switch (am->update) { 4736285c0a3SHansol Suh case TAO_ADMM_UPDATE_BASIC: 4746285c0a3SHansol Suh am->muold = am->mu; 4756285c0a3SHansol Suh break; 4766285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE: 4776285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 4786285c0a3SHansol Suh if (tao->niter == 0) { 4796285c0a3SHansol Suh ierr = VecCopy(am->y, am->y0);CHKERRQ(ierr); 4806285c0a3SHansol Suh ierr = VecWAXPY(am->residual, 1., am->Ax, am->Bzold);CHKERRQ(ierr); 4816285c0a3SHansol Suh if (am->constraint) { 4826285c0a3SHansol Suh ierr = VecAXPY(am->residual, -1., am->constraint);CHKERRQ(ierr); 4836285c0a3SHansol Suh } 4846285c0a3SHansol Suh ierr = VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold);CHKERRQ(ierr); 4856285c0a3SHansol Suh ierr = VecCopy(am->Ax, am->Axold);CHKERRQ(ierr); 4866285c0a3SHansol Suh ierr = VecCopy(am->Bz, am->Bz0);CHKERRQ(ierr); 4876285c0a3SHansol Suh am->muold = am->mu; 4886285c0a3SHansol Suh } else if (tao->niter % am->T == 1) { 4896285c0a3SHansol Suh /* we have compute Bzold in a previous iteration, and we computed Ax above */ 4906285c0a3SHansol Suh ierr = VecWAXPY(am->residual, 1., am->Ax, am->Bzold);CHKERRQ(ierr); 4916285c0a3SHansol Suh if (am->constraint) { 4926285c0a3SHansol Suh ierr = VecAXPY(am->residual, -1., am->constraint);CHKERRQ(ierr); 4936285c0a3SHansol Suh } 4946285c0a3SHansol Suh ierr = VecWAXPY(am->yhat, -am->mu, am->residual, am->yold);CHKERRQ(ierr); 4956285c0a3SHansol Suh ierr = AdaptiveADMMPenaltyUpdate(tao);CHKERRQ(ierr); 4966285c0a3SHansol Suh ierr = VecCopy(am->Ax, am->Axold);CHKERRQ(ierr); 4976285c0a3SHansol Suh ierr = VecCopy(am->Bz, am->Bz0);CHKERRQ(ierr); 4986285c0a3SHansol Suh ierr = VecCopy(am->yhat, am->yhatold);CHKERRQ(ierr); 4996285c0a3SHansol Suh ierr = VecCopy(am->y, am->y0);CHKERRQ(ierr); 5006285c0a3SHansol Suh } else { 5016285c0a3SHansol Suh am->muold = am->mu; 5026285c0a3SHansol Suh } 5036285c0a3SHansol Suh break; 5046285c0a3SHansol Suh default: 5056285c0a3SHansol Suh break; 5066285c0a3SHansol Suh } 5076285c0a3SHansol Suh tao->niter++; 5086285c0a3SHansol Suh 5096285c0a3SHansol Suh /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/ 5106285c0a3SHansol Suh switch (am->regswitch) { 5116285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_USER: 5126285c0a3SHansol Suh if (is_reg_shell) { 5136285c0a3SHansol Suh ierr = ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func);CHKERRQ(ierr); 5146285c0a3SHansol Suh } else { 5156285c0a3SHansol Suh (*am->ops->regobjgrad)(am->subsolverZ,am->subsolverX->solution,®_func,tempL,am->regobjgradP); 5166285c0a3SHansol Suh } 5176285c0a3SHansol Suh break; 5186285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_SOFT_THRESH: 5196285c0a3SHansol Suh ierr = ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func);CHKERRQ(ierr); 5206285c0a3SHansol Suh break; 5216285c0a3SHansol Suh } 5226285c0a3SHansol Suh ierr = VecCopy(am->y,am->yold);CHKERRQ(ierr); 5236285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 5246285c0a3SHansol Suh ierr = VecNorm(am->residual,NORM_2,&am->resnorm);CHKERRQ(ierr); 5256285c0a3SHansol Suh ierr = TaoLogConvergenceHistory(tao,am->last_misfit_val + reg_func,am->dualres,am->resnorm,tao->ksp_its);CHKERRQ(ierr); 5266285c0a3SHansol Suh 5276285c0a3SHansol Suh ierr = TaoMonitor(tao,tao->niter,am->last_misfit_val + reg_func,am->dualres,am->resnorm,1.0);CHKERRQ(ierr); 5286285c0a3SHansol Suh ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 5296285c0a3SHansol Suh } 5306285c0a3SHansol Suh /* Update vectors */ 5316285c0a3SHansol Suh ierr = VecCopy(am->subsolverX->solution,tao->solution);CHKERRQ(ierr); 5326285c0a3SHansol Suh ierr = VecCopy(am->subsolverX->gradient,tao->gradient);CHKERRQ(ierr); 5336285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", NULL);CHKERRQ(ierr); 5346285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", NULL);CHKERRQ(ierr); 5356285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL);CHKERRQ(ierr); 5366285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL);CHKERRQ(ierr); 5376285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL);CHKERRQ(ierr); 5386285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL);CHKERRQ(ierr); 5396285c0a3SHansol Suh PetscFunctionReturn(0); 5406285c0a3SHansol Suh } 5416285c0a3SHansol Suh 5426285c0a3SHansol Suh static PetscErrorCode TaoSetFromOptions_ADMM(PetscOptionItems *PetscOptionsObject,Tao tao) 5436285c0a3SHansol Suh { 5446285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5456285c0a3SHansol Suh PetscErrorCode ierr; 5466285c0a3SHansol Suh 5476285c0a3SHansol Suh PetscFunctionBegin; 5486285c0a3SHansol Suh ierr = PetscOptionsHead(PetscOptionsObject,"ADMM problem that solves f(x) in a form of f(x) + g(z) subject to x - z = 0. Norm 1 and 2 are supported. Different subsolver routines can be selected. ");CHKERRQ(ierr); 5496285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_regularizer_coefficient","regularizer constant","",am->lambda,&am->lambda,NULL);CHKERRQ(ierr); 5506285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_spectral_penalty","Constant for Augmented Lagrangian term.","",am->mu,&am->mu,NULL);CHKERRQ(ierr); 5516285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_relaxation_parameter","x relaxation parameter for Z update.","",am->gamma,&am->gamma,NULL);CHKERRQ(ierr); 5526285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_tolerance_update_factor","ADMM dynamic tolerance update factor.","",am->tol,&am->tol,NULL);CHKERRQ(ierr); 5536285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_spectral_penalty_update_factor","ADMM spectral penalty update curvature safeguard value.","",am->orthval,&am->orthval,NULL);CHKERRQ(ierr); 5546285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_minimum_spectral_penalty","Set ADMM minimum spectral penalty.","",am->mumin,&am->mumin,NULL);CHKERRQ(ierr); 5556285c0a3SHansol Suh ierr = PetscOptionsEnum("-tao_admm_dual_update","Lagrangian dual update policy","TaoADMMUpdateType", 5566285c0a3SHansol Suh TaoADMMUpdateTypes,(PetscEnum)am->update,(PetscEnum*)&am->update,NULL);CHKERRQ(ierr); 5576285c0a3SHansol Suh ierr = PetscOptionsEnum("-tao_admm_regularizer_type","ADMM regularizer update rule","TaoADMMRegularizerType", 5586285c0a3SHansol Suh TaoADMMRegularizerTypes,(PetscEnum)am->regswitch,(PetscEnum*)&am->regswitch,NULL);CHKERRQ(ierr); 5596285c0a3SHansol Suh ierr = PetscOptionsTail();CHKERRQ(ierr); 5606285c0a3SHansol Suh ierr = TaoSetFromOptions(am->subsolverX);CHKERRQ(ierr); 5616285c0a3SHansol Suh if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 5626285c0a3SHansol Suh ierr = TaoSetFromOptions(am->subsolverZ);CHKERRQ(ierr); 5636285c0a3SHansol Suh } 5646285c0a3SHansol Suh PetscFunctionReturn(0); 5656285c0a3SHansol Suh } 5666285c0a3SHansol Suh 5676285c0a3SHansol Suh static PetscErrorCode TaoView_ADMM(Tao tao,PetscViewer viewer) 5686285c0a3SHansol Suh { 5696285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5706285c0a3SHansol Suh PetscErrorCode ierr; 5716285c0a3SHansol Suh 5726285c0a3SHansol Suh PetscFunctionBegin; 5736285c0a3SHansol Suh ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 5746285c0a3SHansol Suh ierr = TaoView(am->subsolverX,viewer);CHKERRQ(ierr); 5756285c0a3SHansol Suh ierr = TaoView(am->subsolverZ,viewer);CHKERRQ(ierr); 5766285c0a3SHansol Suh ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 5776285c0a3SHansol Suh PetscFunctionReturn(0); 5786285c0a3SHansol Suh } 5796285c0a3SHansol Suh 5806285c0a3SHansol Suh static PetscErrorCode TaoSetUp_ADMM(Tao tao) 5816285c0a3SHansol Suh { 5826285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5836285c0a3SHansol Suh PetscErrorCode ierr; 5846285c0a3SHansol Suh PetscInt n,N,M; 5856285c0a3SHansol Suh 5866285c0a3SHansol Suh PetscFunctionBegin; 5876285c0a3SHansol Suh ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 5886285c0a3SHansol Suh ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 5896285c0a3SHansol Suh /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */ 5906285c0a3SHansol Suh if (!am->JB) { 5916285c0a3SHansol Suh am->zJI = PETSC_TRUE; 592ea13f565SStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JB);CHKERRQ(ierr); 5936285c0a3SHansol Suh ierr = MatShellSetOperation(am->JB,MATOP_MULT,(void (*)(void))JacobianIdentityB);CHKERRQ(ierr); 5946285c0a3SHansol Suh ierr = MatShellSetOperation(am->JB,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentityB);CHKERRQ(ierr); 5956285c0a3SHansol Suh ierr = MatShellSetOperation(am->JB,MATOP_TRANSPOSE_MAT_MULT,(void (*)(void))JacobianIdentityB);CHKERRQ(ierr); 5966285c0a3SHansol Suh am->JBpre = am->JB; 5976285c0a3SHansol Suh } 5986285c0a3SHansol Suh if (!am->JA) { 5996285c0a3SHansol Suh am->xJI = PETSC_TRUE; 600ea13f565SStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JA);CHKERRQ(ierr); 6016285c0a3SHansol Suh ierr = MatShellSetOperation(am->JA,MATOP_MULT,(void (*)(void))JacobianIdentity);CHKERRQ(ierr); 6026285c0a3SHansol Suh ierr = MatShellSetOperation(am->JA,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentity);CHKERRQ(ierr); 6036285c0a3SHansol Suh ierr = MatShellSetOperation(am->JA,MATOP_TRANSPOSE_MAT_MULT,(void (*)(void))JacobianIdentity);CHKERRQ(ierr); 6046285c0a3SHansol Suh am->JApre = am->JA; 6056285c0a3SHansol Suh } 6066285c0a3SHansol Suh ierr = MatCreateVecs(am->JA,NULL,&am->Ax);CHKERRQ(ierr); 6076285c0a3SHansol Suh if (!tao->gradient) { 6086285c0a3SHansol Suh ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 6096285c0a3SHansol Suh } 610a82e8c82SStefano Zampini ierr = TaoSetSolution(am->subsolverX, tao->solution);CHKERRQ(ierr); 6116285c0a3SHansol Suh if (!am->z) { 6126285c0a3SHansol Suh ierr = VecDuplicate(tao->solution,&am->z);CHKERRQ(ierr); 6136285c0a3SHansol Suh ierr = VecSet(am->z,0.0);CHKERRQ(ierr); 6146285c0a3SHansol Suh } 615a82e8c82SStefano Zampini ierr = TaoSetSolution(am->subsolverZ, am->z);CHKERRQ(ierr); 6166285c0a3SHansol Suh if (!am->workLeft) { 6176285c0a3SHansol Suh ierr = VecDuplicate(tao->solution,&am->workLeft);CHKERRQ(ierr); 6186285c0a3SHansol Suh } 6196285c0a3SHansol Suh if (!am->Axold) { 6206285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Axold);CHKERRQ(ierr); 6216285c0a3SHansol Suh } 6226285c0a3SHansol Suh if (!am->workJacobianRight) { 6236285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->workJacobianRight);CHKERRQ(ierr); 6246285c0a3SHansol Suh } 6256285c0a3SHansol Suh if (!am->workJacobianRight2) { 6266285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->workJacobianRight2);CHKERRQ(ierr); 6276285c0a3SHansol Suh } 6286285c0a3SHansol Suh if (!am->Bz) { 6296285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Bz);CHKERRQ(ierr); 6306285c0a3SHansol Suh } 6316285c0a3SHansol Suh if (!am->Bzold) { 6326285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Bzold);CHKERRQ(ierr); 6336285c0a3SHansol Suh } 6346285c0a3SHansol Suh if (!am->Bz0) { 6356285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Bz0);CHKERRQ(ierr); 6366285c0a3SHansol Suh } 6376285c0a3SHansol Suh if (!am->y) { 6386285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->y);CHKERRQ(ierr); 6396285c0a3SHansol Suh ierr = VecSet(am->y,0.0);CHKERRQ(ierr); 6406285c0a3SHansol Suh } 6416285c0a3SHansol Suh if (!am->yold) { 6426285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->yold);CHKERRQ(ierr); 6436285c0a3SHansol Suh ierr = VecSet(am->yold,0.0);CHKERRQ(ierr); 6446285c0a3SHansol Suh } 6456285c0a3SHansol Suh if (!am->y0) { 6466285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->y0);CHKERRQ(ierr); 6476285c0a3SHansol Suh ierr = VecSet(am->y0,0.0);CHKERRQ(ierr); 6486285c0a3SHansol Suh } 6496285c0a3SHansol Suh if (!am->yhat) { 6506285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->yhat);CHKERRQ(ierr); 6516285c0a3SHansol Suh ierr = VecSet(am->yhat,0.0);CHKERRQ(ierr); 6526285c0a3SHansol Suh } 6536285c0a3SHansol Suh if (!am->yhatold) { 6546285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->yhatold);CHKERRQ(ierr); 6556285c0a3SHansol Suh ierr = VecSet(am->yhatold,0.0);CHKERRQ(ierr); 6566285c0a3SHansol Suh } 6576285c0a3SHansol Suh if (!am->residual) { 6586285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->residual);CHKERRQ(ierr); 6596285c0a3SHansol Suh ierr = VecSet(am->residual,0.0);CHKERRQ(ierr); 6606285c0a3SHansol Suh } 6616285c0a3SHansol Suh if (!am->constraint) { 6626285c0a3SHansol Suh am->constraint = NULL; 6636285c0a3SHansol Suh } else { 6646285c0a3SHansol Suh ierr = VecGetSize(am->constraint,&M);CHKERRQ(ierr); 665*3c859ba3SBarry Smith PetscCheck(M == N,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Solution vector and constraint vector must be of same size!"); 6666285c0a3SHansol Suh } 6676285c0a3SHansol Suh 6686285c0a3SHansol Suh /* Save changed tao tolerance for adaptive tolerance */ 6696285c0a3SHansol Suh if (tao->gatol_changed) { 6706285c0a3SHansol Suh am->gatol_admm = tao->gatol; 6716285c0a3SHansol Suh } 6726285c0a3SHansol Suh if (tao->catol_changed) { 6736285c0a3SHansol Suh am->catol_admm = tao->catol; 6746285c0a3SHansol Suh } 6756285c0a3SHansol Suh 6766285c0a3SHansol Suh /*Update spectral and dual elements to X subsolver */ 677a82e8c82SStefano Zampini ierr = TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao);CHKERRQ(ierr); 6786285c0a3SHansol Suh ierr = TaoSetJacobianEqualityRoutine(am->subsolverX,am->JA,am->JApre, am->ops->misfitjac, am->misfitjacobianP);CHKERRQ(ierr); 6796285c0a3SHansol Suh ierr = TaoSetJacobianEqualityRoutine(am->subsolverZ,am->JB,am->JBpre, am->ops->regjac, am->regjacobianP);CHKERRQ(ierr); 6806285c0a3SHansol Suh PetscFunctionReturn(0); 6816285c0a3SHansol Suh } 6826285c0a3SHansol Suh 6836285c0a3SHansol Suh static PetscErrorCode TaoDestroy_ADMM(Tao tao) 6846285c0a3SHansol Suh { 6856285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 6866285c0a3SHansol Suh PetscErrorCode ierr; 6876285c0a3SHansol Suh 6886285c0a3SHansol Suh PetscFunctionBegin; 6896285c0a3SHansol Suh ierr = VecDestroy(&am->z);CHKERRQ(ierr); 6906285c0a3SHansol Suh ierr = VecDestroy(&am->Ax);CHKERRQ(ierr); 6916285c0a3SHansol Suh ierr = VecDestroy(&am->Axold);CHKERRQ(ierr); 6926285c0a3SHansol Suh ierr = VecDestroy(&am->Bz);CHKERRQ(ierr); 6936285c0a3SHansol Suh ierr = VecDestroy(&am->Bzold);CHKERRQ(ierr); 6946285c0a3SHansol Suh ierr = VecDestroy(&am->Bz0);CHKERRQ(ierr); 6956285c0a3SHansol Suh ierr = VecDestroy(&am->residual);CHKERRQ(ierr); 6966285c0a3SHansol Suh ierr = VecDestroy(&am->y);CHKERRQ(ierr); 6976285c0a3SHansol Suh ierr = VecDestroy(&am->yold);CHKERRQ(ierr); 6986285c0a3SHansol Suh ierr = VecDestroy(&am->y0);CHKERRQ(ierr); 6996285c0a3SHansol Suh ierr = VecDestroy(&am->yhat);CHKERRQ(ierr); 7006285c0a3SHansol Suh ierr = VecDestroy(&am->yhatold);CHKERRQ(ierr); 7016285c0a3SHansol Suh ierr = VecDestroy(&am->workLeft);CHKERRQ(ierr); 7026285c0a3SHansol Suh ierr = VecDestroy(&am->workJacobianRight);CHKERRQ(ierr); 7036285c0a3SHansol Suh ierr = VecDestroy(&am->workJacobianRight2);CHKERRQ(ierr); 7046285c0a3SHansol Suh 7056285c0a3SHansol Suh ierr = MatDestroy(&am->JA);CHKERRQ(ierr); 7066285c0a3SHansol Suh ierr = MatDestroy(&am->JB);CHKERRQ(ierr); 7076285c0a3SHansol Suh if (!am->xJI) { 7086285c0a3SHansol Suh ierr = MatDestroy(&am->JApre);CHKERRQ(ierr); 7096285c0a3SHansol Suh } 7106285c0a3SHansol Suh if (!am->zJI) { 7116285c0a3SHansol Suh ierr = MatDestroy(&am->JBpre);CHKERRQ(ierr); 7126285c0a3SHansol Suh } 7136285c0a3SHansol Suh if (am->Hx) { 7146285c0a3SHansol Suh ierr = MatDestroy(&am->Hx);CHKERRQ(ierr); 7156285c0a3SHansol Suh ierr = MatDestroy(&am->Hxpre);CHKERRQ(ierr); 7166285c0a3SHansol Suh } 7176285c0a3SHansol Suh if (am->Hz) { 7186285c0a3SHansol Suh ierr = MatDestroy(&am->Hz);CHKERRQ(ierr); 7196285c0a3SHansol Suh ierr = MatDestroy(&am->Hzpre);CHKERRQ(ierr); 7206285c0a3SHansol Suh } 7216285c0a3SHansol Suh ierr = MatDestroy(&am->ATA);CHKERRQ(ierr); 7226285c0a3SHansol Suh ierr = MatDestroy(&am->BTB);CHKERRQ(ierr); 7236285c0a3SHansol Suh ierr = TaoDestroy(&am->subsolverX);CHKERRQ(ierr); 7246285c0a3SHansol Suh ierr = TaoDestroy(&am->subsolverZ);CHKERRQ(ierr); 7256285c0a3SHansol Suh am->parent = NULL; 7266285c0a3SHansol Suh ierr = PetscFree(tao->data);CHKERRQ(ierr); 7276285c0a3SHansol Suh PetscFunctionReturn(0); 7286285c0a3SHansol Suh } 7296285c0a3SHansol Suh 7306285c0a3SHansol Suh /*MC 7316285c0a3SHansol Suh 7326285c0a3SHansol Suh TAOADMM - Alternating direction method of multipliers method fo solving linear problems with 7336285c0a3SHansol Suh constraints. in a min_x f(x) + g(z) s.t. Ax+Bz=c. 734a5b23f4aSJose E. Roman This algorithm employs two sub Tao solvers, of which type can be specified 7356285c0a3SHansol Suh by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers. 7366285c0a3SHansol Suh Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via 7376285c0a3SHansol Suh TaoADMMSet{Misfit,Regularizer}HessianChangeStatus. 7386285c0a3SHansol Suh Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type. 7396285c0a3SHansol Suh There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update, 740a5b23f4aSJose E. Roman currently there are basic option and adaptive option. 7416285c0a3SHansol Suh Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian. 7426285c0a3SHansol Suh c can be set with TaoADMMSetConstraintVectorRHS. 7436285c0a3SHansol Suh The user can also provide regularizer weight for second subsolver. 7446285c0a3SHansol Suh 7456285c0a3SHansol Suh References: 746606c0280SSatish Balay . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom 7476285c0a3SHansol Suh "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation" 7486285c0a3SHansol Suh The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017. 7496285c0a3SHansol Suh 7506285c0a3SHansol Suh Options Database Keys: 7516285c0a3SHansol Suh + -tao_admm_regularizer_coefficient - regularizer constant (default 1.e-6) 7526285c0a3SHansol Suh . -tao_admm_spectral_penalty - Constant for Augmented Lagrangian term (default 1.) 7536285c0a3SHansol Suh . -tao_admm_relaxation_parameter - relaxation parameter for Z update (default 1.) 7546285c0a3SHansol Suh . -tao_admm_tolerance_update_factor - ADMM dynamic tolerance update factor (default 1.e-12) 7556285c0a3SHansol Suh . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2) 7566285c0a3SHansol Suh . -tao_admm_minimum_spectral_penalty - Set ADMM minimum spectral penalty (default 0) 7576285c0a3SHansol Suh . -tao_admm_dual_update - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic") 7586285c0a3SHansol Suh - -tao_admm_regularizer_type - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold") 7596285c0a3SHansol Suh 7606285c0a3SHansol Suh Level: beginner 7616285c0a3SHansol Suh 7626285c0a3SHansol Suh .seealso: TaoADMMSetMisfitHessianChangeStatus(), TaoADMMSetRegHessianChangeStatus(), TaoADMMGetSpectralPenalty(), 7636285c0a3SHansol Suh TaoADMMGetMisfitSubsolver(), TaoADMMGetRegularizationSubsolver(), TaoADMMSetConstraintVectorRHS(), 7646285c0a3SHansol Suh TaoADMMSetMinimumSpectralPenalty(), TaoADMMSetRegularizerCoefficient(), 7656285c0a3SHansol Suh TaoADMMSetRegularizerConstraintJacobian(), TaoADMMSetMisfitConstraintJacobian(), 7666285c0a3SHansol Suh TaoADMMSetMisfitObjectiveAndGradientRoutine(), TaoADMMSetMisfitHessianRoutine(), 7676285c0a3SHansol Suh TaoADMMSetRegularizerObjectiveAndGradientRoutine(), TaoADMMSetRegularizerHessianRoutine(), 7686285c0a3SHansol Suh TaoGetADMMParentTao(), TaoADMMGetDualVector(), TaoADMMSetRegularizerType(), 7696285c0a3SHansol Suh TaoADMMGetRegularizerType(), TaoADMMSetUpdateType(), TaoADMMGetUpdateType() 7706285c0a3SHansol Suh M*/ 7716285c0a3SHansol Suh 7726285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao) 7736285c0a3SHansol Suh { 7746285c0a3SHansol Suh TAO_ADMM *am; 7756285c0a3SHansol Suh PetscErrorCode ierr; 7766285c0a3SHansol Suh 7776285c0a3SHansol Suh PetscFunctionBegin; 7786285c0a3SHansol Suh ierr = PetscNewLog(tao,&am);CHKERRQ(ierr); 7796285c0a3SHansol Suh 7806285c0a3SHansol Suh tao->ops->destroy = TaoDestroy_ADMM; 7816285c0a3SHansol Suh tao->ops->setup = TaoSetUp_ADMM; 7826285c0a3SHansol Suh tao->ops->setfromoptions = TaoSetFromOptions_ADMM; 7836285c0a3SHansol Suh tao->ops->view = TaoView_ADMM; 7846285c0a3SHansol Suh tao->ops->solve = TaoSolve_ADMM; 7856285c0a3SHansol Suh 7866285c0a3SHansol Suh tao->data = (void*)am; 7876285c0a3SHansol Suh am->l1epsilon = 1e-6; 7886285c0a3SHansol Suh am->lambda = 1e-4; 7896285c0a3SHansol Suh am->mu = 1.; 7906285c0a3SHansol Suh am->muold = 0.; 7916285c0a3SHansol Suh am->mueps = PETSC_MACHINE_EPSILON; 7926285c0a3SHansol Suh am->mumin = 0.; 7936285c0a3SHansol Suh am->orthval = 0.2; 7946285c0a3SHansol Suh am->T = 2; 7956285c0a3SHansol Suh am->parent = tao; 7966285c0a3SHansol Suh am->update = TAO_ADMM_UPDATE_BASIC; 7976285c0a3SHansol Suh am->regswitch = TAO_ADMM_REGULARIZER_SOFT_THRESH; 7986285c0a3SHansol Suh am->tol = PETSC_SMALL; 7996285c0a3SHansol Suh am->const_norm = 0; 8006285c0a3SHansol Suh am->resnorm = 0; 8016285c0a3SHansol Suh am->dualres = 0; 80283c8fe1dSLisandro Dalcin am->ops->regobjgrad = NULL; 80383c8fe1dSLisandro Dalcin am->ops->reghess = NULL; 8046285c0a3SHansol Suh am->gamma = 1; 80583c8fe1dSLisandro Dalcin am->regobjgradP = NULL; 80683c8fe1dSLisandro Dalcin am->reghessP = NULL; 8076285c0a3SHansol Suh am->gatol_admm = 1e-8; 8086285c0a3SHansol Suh am->catol_admm = 0; 8096285c0a3SHansol Suh am->Hxchange = PETSC_TRUE; 8106285c0a3SHansol Suh am->Hzchange = PETSC_TRUE; 8116285c0a3SHansol Suh am->Hzbool = PETSC_TRUE; 8126285c0a3SHansol Suh am->Hxbool = PETSC_TRUE; 8136285c0a3SHansol Suh 8146285c0a3SHansol Suh ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverX);CHKERRQ(ierr); 8156285c0a3SHansol Suh ierr = TaoSetOptionsPrefix(am->subsolverX,"misfit_");CHKERRQ(ierr); 8166285c0a3SHansol Suh ierr = PetscObjectIncrementTabLevel((PetscObject)am->subsolverX,(PetscObject)tao,1);CHKERRQ(ierr); 8176285c0a3SHansol Suh ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverZ);CHKERRQ(ierr); 8186285c0a3SHansol Suh ierr = TaoSetOptionsPrefix(am->subsolverZ,"reg_");CHKERRQ(ierr); 8196285c0a3SHansol Suh ierr = PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ,(PetscObject)tao,1);CHKERRQ(ierr); 8206285c0a3SHansol Suh 8216285c0a3SHansol Suh ierr = TaoSetType(am->subsolverX,TAONLS);CHKERRQ(ierr); 8226285c0a3SHansol Suh ierr = TaoSetType(am->subsolverZ,TAONLS);CHKERRQ(ierr); 8236285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", (PetscObject) tao);CHKERRQ(ierr); 8246285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", (PetscObject) tao);CHKERRQ(ierr); 8256285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",TaoADMMSetRegularizerType_ADMM);CHKERRQ(ierr); 8266285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",TaoADMMGetRegularizerType_ADMM);CHKERRQ(ierr); 8276285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",TaoADMMSetUpdateType_ADMM);CHKERRQ(ierr); 8286285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",TaoADMMGetUpdateType_ADMM);CHKERRQ(ierr); 8296285c0a3SHansol Suh PetscFunctionReturn(0); 8306285c0a3SHansol Suh } 8316285c0a3SHansol Suh 8326285c0a3SHansol Suh /*@ 8336285c0a3SHansol Suh TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines whether Hessian matrix of misfit subsolver changes with respect to input vector. 8346285c0a3SHansol Suh 8356285c0a3SHansol Suh Collective on Tao 8366285c0a3SHansol Suh 8376285c0a3SHansol Suh Input Parameters: 8387f5c9be9SBarry Smith + tao - the Tao solver context. 8397f5c9be9SBarry Smith - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 8406285c0a3SHansol Suh 8416285c0a3SHansol Suh Level: advanced 842a5a2f7acSBarry Smith 843a5a2f7acSBarry Smith .seealso: TAOADMM 844a5a2f7acSBarry Smith 8456285c0a3SHansol Suh @*/ 8466285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b) 8476285c0a3SHansol Suh { 8486285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8496285c0a3SHansol Suh 8506285c0a3SHansol Suh PetscFunctionBegin; 8516285c0a3SHansol Suh am->Hxchange = b; 8526285c0a3SHansol Suh PetscFunctionReturn(0); 8536285c0a3SHansol Suh } 8546285c0a3SHansol Suh 8556285c0a3SHansol Suh /*@ 8566285c0a3SHansol Suh TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector. 8577f5c9be9SBarry Smith 8586285c0a3SHansol Suh Collective on Tao 8596285c0a3SHansol Suh 8606285c0a3SHansol Suh Input Parameters: 8616285c0a3SHansol Suh + tao - the Tao solver context 8627f5c9be9SBarry Smith - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 8636285c0a3SHansol Suh 8646285c0a3SHansol Suh Level: advanced 865a5a2f7acSBarry Smith 866a5a2f7acSBarry Smith .seealso: TAOADMM 867a5a2f7acSBarry Smith 8686285c0a3SHansol Suh @*/ 8696285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b) 8706285c0a3SHansol Suh { 8716285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8726285c0a3SHansol Suh 8736285c0a3SHansol Suh PetscFunctionBegin; 8746285c0a3SHansol Suh am->Hzchange = b; 8756285c0a3SHansol Suh PetscFunctionReturn(0); 8766285c0a3SHansol Suh } 8776285c0a3SHansol Suh 8786285c0a3SHansol Suh /*@ 8796285c0a3SHansol Suh TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value 8806285c0a3SHansol Suh 8816285c0a3SHansol Suh Collective on Tao 8826285c0a3SHansol Suh 8836285c0a3SHansol Suh Input Parameters: 8846285c0a3SHansol Suh + tao - the Tao solver context 8857f5c9be9SBarry Smith - mu - spectral penalty 8866285c0a3SHansol Suh 8876285c0a3SHansol Suh Level: advanced 8886285c0a3SHansol Suh 889a5a2f7acSBarry Smith .seealso: TaoADMMSetMinimumSpectralPenalty(), TAOADMM 8906285c0a3SHansol Suh @*/ 8916285c0a3SHansol Suh PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu) 8926285c0a3SHansol Suh { 8936285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8946285c0a3SHansol Suh 8956285c0a3SHansol Suh PetscFunctionBegin; 8966285c0a3SHansol Suh am->mu = mu; 8976285c0a3SHansol Suh PetscFunctionReturn(0); 8986285c0a3SHansol Suh } 8996285c0a3SHansol Suh 9006285c0a3SHansol Suh /*@ 9016285c0a3SHansol Suh TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value 9026285c0a3SHansol Suh 9036285c0a3SHansol Suh Collective on Tao 9046285c0a3SHansol Suh 9057f5c9be9SBarry Smith Input Parameter: 9067f5c9be9SBarry Smith . tao - the Tao solver context 9077f5c9be9SBarry Smith 9087f5c9be9SBarry Smith Output Parameter: 9097f5c9be9SBarry Smith . mu - spectral penalty 9106285c0a3SHansol Suh 9116285c0a3SHansol Suh Level: advanced 9126285c0a3SHansol Suh 913a5a2f7acSBarry Smith .seealso: TaoADMMSetMinimumSpectralPenalty(), TaoADMMSetSpectralPenalty(), TAOADMM 9146285c0a3SHansol Suh @*/ 9156285c0a3SHansol Suh PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu) 9166285c0a3SHansol Suh { 9176285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9186285c0a3SHansol Suh 9196285c0a3SHansol Suh PetscFunctionBegin; 9206285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 9216285c0a3SHansol Suh PetscValidRealPointer(mu,2); 9226285c0a3SHansol Suh *mu = am->mu; 9236285c0a3SHansol Suh PetscFunctionReturn(0); 9246285c0a3SHansol Suh } 9256285c0a3SHansol Suh 9266285c0a3SHansol Suh /*@ 9276285c0a3SHansol Suh TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM 9286285c0a3SHansol Suh 9296285c0a3SHansol Suh Collective on Tao 9306285c0a3SHansol Suh 9317f5c9be9SBarry Smith Input Parameter: 9327f5c9be9SBarry Smith . tao - the Tao solver context 9337f5c9be9SBarry Smith 9347f5c9be9SBarry Smith Output Parameter: 9357f5c9be9SBarry Smith . misfit - the Tao subsolver context 9366285c0a3SHansol Suh 9376285c0a3SHansol Suh Level: advanced 938a5a2f7acSBarry Smith 939a5a2f7acSBarry Smith .seealso: TAOADMM 940a5a2f7acSBarry Smith 9416285c0a3SHansol Suh @*/ 9426285c0a3SHansol Suh PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit) 9436285c0a3SHansol Suh { 9446285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9456285c0a3SHansol Suh 9466285c0a3SHansol Suh PetscFunctionBegin; 9476285c0a3SHansol Suh *misfit = am->subsolverX; 9486285c0a3SHansol Suh PetscFunctionReturn(0); 9496285c0a3SHansol Suh } 9506285c0a3SHansol Suh 9516285c0a3SHansol Suh /*@ 9526285c0a3SHansol Suh TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM 9536285c0a3SHansol Suh 9546285c0a3SHansol Suh Collective on Tao 9556285c0a3SHansol Suh 9567f5c9be9SBarry Smith Input Parameter: 9577f5c9be9SBarry Smith . tao - the Tao solver context 9587f5c9be9SBarry Smith 9597f5c9be9SBarry Smith Output Parameter: 9607f5c9be9SBarry Smith . reg - the Tao subsolver context 9616285c0a3SHansol Suh 9626285c0a3SHansol Suh Level: advanced 963a5a2f7acSBarry Smith 964a5a2f7acSBarry Smith .seealso: TAOADMM 965a5a2f7acSBarry Smith 9666285c0a3SHansol Suh @*/ 9676285c0a3SHansol Suh PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg) 9686285c0a3SHansol Suh { 9696285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9706285c0a3SHansol Suh 9716285c0a3SHansol Suh PetscFunctionBegin; 9726285c0a3SHansol Suh *reg = am->subsolverZ; 9736285c0a3SHansol Suh PetscFunctionReturn(0); 9746285c0a3SHansol Suh } 9756285c0a3SHansol Suh 9766285c0a3SHansol Suh /*@ 9776285c0a3SHansol Suh TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM 9786285c0a3SHansol Suh 9796285c0a3SHansol Suh Collective on Tao 9806285c0a3SHansol Suh 9816285c0a3SHansol Suh Input Parameters: 9826285c0a3SHansol Suh + tao - the Tao solver context 9836285c0a3SHansol Suh - c - RHS vector 9846285c0a3SHansol Suh 9856285c0a3SHansol Suh Level: advanced 986a5a2f7acSBarry Smith 987a5a2f7acSBarry Smith .seealso: TAOADMM 988a5a2f7acSBarry Smith 9896285c0a3SHansol Suh @*/ 9906285c0a3SHansol Suh PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c) 9916285c0a3SHansol Suh { 9926285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9936285c0a3SHansol Suh 9946285c0a3SHansol Suh PetscFunctionBegin; 9956285c0a3SHansol Suh am->constraint = c; 9966285c0a3SHansol Suh PetscFunctionReturn(0); 9976285c0a3SHansol Suh } 9986285c0a3SHansol Suh 9996285c0a3SHansol Suh /*@ 10007f5c9be9SBarry Smith TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty 10016285c0a3SHansol Suh 10026285c0a3SHansol Suh Collective on Tao 10036285c0a3SHansol Suh 10046285c0a3SHansol Suh Input Parameters: 10056285c0a3SHansol Suh + tao - the Tao solver context 10066285c0a3SHansol Suh - mu - minimum spectral penalty value 10076285c0a3SHansol Suh 10086285c0a3SHansol Suh Level: advanced 10096285c0a3SHansol Suh 1010a5a2f7acSBarry Smith .seealso: TaoADMMGetSpectralPenalty(), TAOADMM 10116285c0a3SHansol Suh @*/ 10126285c0a3SHansol Suh PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu) 10136285c0a3SHansol Suh { 10146285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10156285c0a3SHansol Suh 10166285c0a3SHansol Suh PetscFunctionBegin; 10176285c0a3SHansol Suh am->mumin= mu; 10186285c0a3SHansol Suh PetscFunctionReturn(0); 10196285c0a3SHansol Suh } 10206285c0a3SHansol Suh 10216285c0a3SHansol Suh /*@ 10226285c0a3SHansol Suh TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case 10236285c0a3SHansol Suh 10246285c0a3SHansol Suh Collective on Tao 10256285c0a3SHansol Suh 10266285c0a3SHansol Suh Input Parameters: 10276285c0a3SHansol Suh + tao - the Tao solver context 10286285c0a3SHansol Suh - lambda - L1-norm regularizer coefficient 10296285c0a3SHansol Suh 10306285c0a3SHansol Suh Level: advanced 10317f5c9be9SBarry Smith 1032a5a2f7acSBarry Smith .seealso: TaoADMMSetMisfitConstraintJacobian(), TaoADMMSetRegularizerConstraintJacobian(), TAOADMM 10337f5c9be9SBarry Smith 10346285c0a3SHansol Suh @*/ 10356285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda) 10366285c0a3SHansol Suh { 10376285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10386285c0a3SHansol Suh 10396285c0a3SHansol Suh PetscFunctionBegin; 10406285c0a3SHansol Suh am->lambda = lambda; 10416285c0a3SHansol Suh PetscFunctionReturn(0); 10426285c0a3SHansol Suh } 10436285c0a3SHansol Suh 10446285c0a3SHansol Suh /*@C 10457f5c9be9SBarry Smith TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable. 10466285c0a3SHansol Suh 10476285c0a3SHansol Suh Collective on Tao 10486285c0a3SHansol Suh 1049a5a2f7acSBarry Smith Input Parameters: 1050a5a2f7acSBarry Smith + tao - the Tao solver context 10516285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix 10526285c0a3SHansol Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 10536285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function 10546285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 10556285c0a3SHansol Suh 10566285c0a3SHansol Suh Level: advanced 10577f5c9be9SBarry Smith 1058a5a2f7acSBarry Smith .seealso: TaoADMMSetRegularizerCoefficient(), TaoADMMSetRegularizerConstraintJacobian(), TAOADMM 10597f5c9be9SBarry Smith 10606285c0a3SHansol Suh @*/ 10616285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 10626285c0a3SHansol Suh { 10636285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10646285c0a3SHansol Suh PetscErrorCode ierr; 10656285c0a3SHansol Suh 10666285c0a3SHansol Suh PetscFunctionBegin; 10676285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 10686285c0a3SHansol Suh if (J) { 10696285c0a3SHansol Suh PetscValidHeaderSpecific(J,MAT_CLASSID,2); 10706285c0a3SHansol Suh PetscCheckSameComm(tao,1,J,2); 10716285c0a3SHansol Suh } 10726285c0a3SHansol Suh if (Jpre) { 10736285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 10746285c0a3SHansol Suh PetscCheckSameComm(tao,1,Jpre,3); 10756285c0a3SHansol Suh } 10766285c0a3SHansol Suh if (ctx) am->misfitjacobianP = ctx; 10776285c0a3SHansol Suh if (func) am->ops->misfitjac = func; 10786285c0a3SHansol Suh 10796285c0a3SHansol Suh if (J) { 10806285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 10816285c0a3SHansol Suh ierr = MatDestroy(&am->JA);CHKERRQ(ierr); 10826285c0a3SHansol Suh am->JA = J; 10836285c0a3SHansol Suh } 10846285c0a3SHansol Suh if (Jpre) { 10856285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 10866285c0a3SHansol Suh ierr = MatDestroy(&am->JApre);CHKERRQ(ierr); 10876285c0a3SHansol Suh am->JApre = Jpre; 10886285c0a3SHansol Suh } 10896285c0a3SHansol Suh PetscFunctionReturn(0); 10906285c0a3SHansol Suh } 10916285c0a3SHansol Suh 10926285c0a3SHansol Suh /*@C 10936285c0a3SHansol Suh TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable. 10946285c0a3SHansol Suh 10956285c0a3SHansol Suh Collective on Tao 10966285c0a3SHansol Suh 10976285c0a3SHansol Suh Input Parameters: 10986285c0a3SHansol Suh + tao - the Tao solver context 10996285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix 11006285c0a3SHansol Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 11016285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function 11026285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 11036285c0a3SHansol Suh 11046285c0a3SHansol Suh Level: advanced 11057f5c9be9SBarry Smith 1106a5a2f7acSBarry Smith .seealso: TaoADMMSetRegularizerCoefficient(), TaoADMMSetMisfitConstraintJacobian(), TAOADMM 11077f5c9be9SBarry Smith 11086285c0a3SHansol Suh @*/ 11096285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 11106285c0a3SHansol Suh { 11116285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11126285c0a3SHansol Suh PetscErrorCode ierr; 11136285c0a3SHansol Suh 11146285c0a3SHansol Suh PetscFunctionBegin; 11156285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11166285c0a3SHansol Suh if (J) { 11176285c0a3SHansol Suh PetscValidHeaderSpecific(J,MAT_CLASSID,2); 11186285c0a3SHansol Suh PetscCheckSameComm(tao,1,J,2); 11196285c0a3SHansol Suh } 11206285c0a3SHansol Suh if (Jpre) { 11216285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 11226285c0a3SHansol Suh PetscCheckSameComm(tao,1,Jpre,3); 11236285c0a3SHansol Suh } 11246285c0a3SHansol Suh if (ctx) am->regjacobianP = ctx; 11256285c0a3SHansol Suh if (func) am->ops->regjac = func; 11266285c0a3SHansol Suh 11276285c0a3SHansol Suh if (J) { 11286285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 11296285c0a3SHansol Suh ierr = MatDestroy(&am->JB);CHKERRQ(ierr); 11306285c0a3SHansol Suh am->JB = J; 11316285c0a3SHansol Suh } 11326285c0a3SHansol Suh if (Jpre) { 11336285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 11346285c0a3SHansol Suh ierr = MatDestroy(&am->JBpre);CHKERRQ(ierr); 11356285c0a3SHansol Suh am->JBpre = Jpre; 11366285c0a3SHansol Suh } 11376285c0a3SHansol Suh PetscFunctionReturn(0); 11386285c0a3SHansol Suh } 11396285c0a3SHansol Suh 11406285c0a3SHansol Suh /*@C 11417f5c9be9SBarry Smith TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function 11427f5c9be9SBarry Smith 11437f5c9be9SBarry Smith Collective on tao 11446285c0a3SHansol Suh 11456285c0a3SHansol Suh Input Parameters: 11466285c0a3SHansol Suh + tao - the Tao context 11476285c0a3SHansol Suh . func - function pointer for the misfit value and gradient evaluation 11486285c0a3SHansol Suh - ctx - user context for the misfit 11496285c0a3SHansol Suh 11506285c0a3SHansol Suh Level: advanced 11517f5c9be9SBarry Smith 1152a5a2f7acSBarry Smith .seealso: TAOADMM 11537f5c9be9SBarry Smith 11546285c0a3SHansol Suh @*/ 11556285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 11566285c0a3SHansol Suh { 11576285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11586285c0a3SHansol Suh 11596285c0a3SHansol Suh PetscFunctionBegin; 11606285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11616285c0a3SHansol Suh am->misfitobjgradP = ctx; 11626285c0a3SHansol Suh am->ops->misfitobjgrad = func; 11636285c0a3SHansol Suh PetscFunctionReturn(0); 11646285c0a3SHansol Suh } 11656285c0a3SHansol Suh 11666285c0a3SHansol Suh /*@C 11676285c0a3SHansol Suh TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back 11686285c0a3SHansol Suh function into the algorithm, to be used for subsolverX. 11696285c0a3SHansol Suh 11707f5c9be9SBarry Smith Collective on tao 11717f5c9be9SBarry Smith 11726285c0a3SHansol Suh Input Parameters: 11736285c0a3SHansol Suh + tao - the Tao context 11746285c0a3SHansol Suh . H - user-created matrix for the Hessian of the misfit term 11756285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term 11766285c0a3SHansol Suh . func - function pointer for the misfit Hessian evaluation 11776285c0a3SHansol Suh - ctx - user context for the misfit Hessian 11786285c0a3SHansol Suh 11796285c0a3SHansol Suh Level: advanced 1180a5a2f7acSBarry Smith 1181a5a2f7acSBarry Smith .seealso: TAOADMM 1182a5a2f7acSBarry Smith 11836285c0a3SHansol Suh @*/ 11846285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 11856285c0a3SHansol Suh { 11866285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11876285c0a3SHansol Suh PetscErrorCode ierr; 11886285c0a3SHansol Suh 11896285c0a3SHansol Suh PetscFunctionBegin; 11906285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11916285c0a3SHansol Suh if (H) { 11926285c0a3SHansol Suh PetscValidHeaderSpecific(H,MAT_CLASSID,2); 11936285c0a3SHansol Suh PetscCheckSameComm(tao,1,H,2); 11946285c0a3SHansol Suh } 11956285c0a3SHansol Suh if (Hpre) { 11966285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 11976285c0a3SHansol Suh PetscCheckSameComm(tao,1,Hpre,3); 11986285c0a3SHansol Suh } 11996285c0a3SHansol Suh if (ctx) { 12006285c0a3SHansol Suh am->misfithessP = ctx; 12016285c0a3SHansol Suh } 12026285c0a3SHansol Suh if (func) { 12036285c0a3SHansol Suh am->ops->misfithess = func; 12046285c0a3SHansol Suh } 12056285c0a3SHansol Suh if (H) { 12066285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr); 12076285c0a3SHansol Suh ierr = MatDestroy(&am->Hx);CHKERRQ(ierr); 12086285c0a3SHansol Suh am->Hx = H; 12096285c0a3SHansol Suh } 12106285c0a3SHansol Suh if (Hpre) { 12116285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr); 12126285c0a3SHansol Suh ierr = MatDestroy(&am->Hxpre);CHKERRQ(ierr); 12136285c0a3SHansol Suh am->Hxpre = Hpre; 12146285c0a3SHansol Suh } 12156285c0a3SHansol Suh PetscFunctionReturn(0); 12166285c0a3SHansol Suh } 12176285c0a3SHansol Suh 12186285c0a3SHansol Suh /*@C 12197f5c9be9SBarry Smith TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function 12207f5c9be9SBarry Smith 12217f5c9be9SBarry Smith Collective on tao 12226285c0a3SHansol Suh 12236285c0a3SHansol Suh Input Parameters: 12246285c0a3SHansol Suh + tao - the Tao context 12256285c0a3SHansol Suh . func - function pointer for the regularizer value and gradient evaluation 12266285c0a3SHansol Suh - ctx - user context for the regularizer 12276285c0a3SHansol Suh 12286285c0a3SHansol Suh Level: advanced 1229a5a2f7acSBarry Smith 1230a5a2f7acSBarry Smith .seealso: TAOADMM 1231a5a2f7acSBarry Smith 12326285c0a3SHansol Suh @*/ 12336285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 12346285c0a3SHansol Suh { 12356285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 12366285c0a3SHansol Suh 12376285c0a3SHansol Suh PetscFunctionBegin; 12386285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 12396285c0a3SHansol Suh am->regobjgradP = ctx; 12406285c0a3SHansol Suh am->ops->regobjgrad = func; 12416285c0a3SHansol Suh PetscFunctionReturn(0); 12426285c0a3SHansol Suh } 12436285c0a3SHansol Suh 12446285c0a3SHansol Suh /*@C 12456285c0a3SHansol Suh TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back 12467f5c9be9SBarry Smith function, to be used for subsolverZ. 12477f5c9be9SBarry Smith 12487f5c9be9SBarry Smith Collective on tao 12496285c0a3SHansol Suh 12506285c0a3SHansol Suh Input Parameters: 12516285c0a3SHansol Suh + tao - the Tao context 12526285c0a3SHansol Suh . H - user-created matrix for the Hessian of the regularization term 12536285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term 12546285c0a3SHansol Suh . func - function pointer for the regularizer Hessian evaluation 12556285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 12566285c0a3SHansol Suh 12576285c0a3SHansol Suh Level: advanced 1258a5a2f7acSBarry Smith 1259a5a2f7acSBarry Smith .seealso: TAOADMM 1260a5a2f7acSBarry Smith 12616285c0a3SHansol Suh @*/ 12626285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 12636285c0a3SHansol Suh { 12646285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 12656285c0a3SHansol Suh PetscErrorCode ierr; 12666285c0a3SHansol Suh 12676285c0a3SHansol Suh PetscFunctionBegin; 12686285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 12696285c0a3SHansol Suh if (H) { 12706285c0a3SHansol Suh PetscValidHeaderSpecific(H,MAT_CLASSID,2); 12716285c0a3SHansol Suh PetscCheckSameComm(tao,1,H,2); 12726285c0a3SHansol Suh } 12736285c0a3SHansol Suh if (Hpre) { 12746285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 12756285c0a3SHansol Suh PetscCheckSameComm(tao,1,Hpre,3); 12766285c0a3SHansol Suh } 12776285c0a3SHansol Suh if (ctx) { 12786285c0a3SHansol Suh am->reghessP = ctx; 12796285c0a3SHansol Suh } 12806285c0a3SHansol Suh if (func) { 12816285c0a3SHansol Suh am->ops->reghess = func; 12826285c0a3SHansol Suh } 12836285c0a3SHansol Suh if (H) { 12846285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr); 12856285c0a3SHansol Suh ierr = MatDestroy(&am->Hz);CHKERRQ(ierr); 12866285c0a3SHansol Suh am->Hz = H; 12876285c0a3SHansol Suh } 12886285c0a3SHansol Suh if (Hpre) { 12896285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr); 12906285c0a3SHansol Suh ierr = MatDestroy(&am->Hzpre);CHKERRQ(ierr); 12916285c0a3SHansol Suh am->Hzpre = Hpre; 12926285c0a3SHansol Suh } 12936285c0a3SHansol Suh PetscFunctionReturn(0); 12946285c0a3SHansol Suh } 12956285c0a3SHansol Suh 12966285c0a3SHansol Suh /*@ 12977f5c9be9SBarry Smith TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver. 12986285c0a3SHansol Suh 12997f5c9be9SBarry Smith Collective on tao 13007f5c9be9SBarry Smith 13017f5c9be9SBarry Smith Input Parameter: 13027f5c9be9SBarry Smith . tao - the Tao context 13037f5c9be9SBarry Smith 13047f5c9be9SBarry Smith Output Parameter: 13057f5c9be9SBarry Smith . admm_tao - the parent Tao context 13066285c0a3SHansol Suh 13076285c0a3SHansol Suh Level: advanced 1308a5a2f7acSBarry Smith 1309a5a2f7acSBarry Smith .seealso: TAOADMM 1310a5a2f7acSBarry Smith 13116285c0a3SHansol Suh @*/ 13126285c0a3SHansol Suh PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao) 13136285c0a3SHansol Suh { 13146285c0a3SHansol Suh PetscErrorCode ierr; 13156285c0a3SHansol Suh 13166285c0a3SHansol Suh PetscFunctionBegin; 13176285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13186285c0a3SHansol Suh ierr = PetscObjectQuery((PetscObject)tao,"TaoGetADMMParentTao_ADMM", (PetscObject*) admm_tao);CHKERRQ(ierr); 13196285c0a3SHansol Suh PetscFunctionReturn(0); 13206285c0a3SHansol Suh } 13216285c0a3SHansol Suh 13226285c0a3SHansol Suh /*@ 13237f5c9be9SBarry Smith TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state 13246285c0a3SHansol Suh 13256285c0a3SHansol Suh Not Collective 13266285c0a3SHansol Suh 13276285c0a3SHansol Suh Input Parameter: 13287f5c9be9SBarry Smith . tao - the Tao context 13297f5c9be9SBarry Smith 13307f5c9be9SBarry Smith Output Parameter: 13317f5c9be9SBarry Smith . Y - the current solution 13326285c0a3SHansol Suh 13336285c0a3SHansol Suh Level: intermediate 1334a5a2f7acSBarry Smith 1335a5a2f7acSBarry Smith .seealso: TAOADMM 1336a5a2f7acSBarry Smith 13376285c0a3SHansol Suh @*/ 13386285c0a3SHansol Suh PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y) 13396285c0a3SHansol Suh { 13406285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 13416285c0a3SHansol Suh 13426285c0a3SHansol Suh PetscFunctionBegin; 13436285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13446285c0a3SHansol Suh *Y = am->y; 13456285c0a3SHansol Suh PetscFunctionReturn(0); 13466285c0a3SHansol Suh } 13476285c0a3SHansol Suh 13486285c0a3SHansol Suh /*@ 13496285c0a3SHansol Suh TaoADMMSetRegularizerType - Set regularizer type for ADMM routine 13506285c0a3SHansol Suh 13516285c0a3SHansol Suh Not Collective 13526285c0a3SHansol Suh 1353d8d19677SJose E. Roman Input Parameters: 13546285c0a3SHansol Suh + tao - the Tao context 13557f5c9be9SBarry Smith - type - regularizer type 13567f5c9be9SBarry Smith 13577f5c9be9SBarry Smith Options Database: 1358147403d9SBarry Smith . -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer 13596285c0a3SHansol Suh 13606285c0a3SHansol Suh Level: intermediate 13616285c0a3SHansol Suh 1362a5a2f7acSBarry Smith .seealso: TaoADMMGetRegularizerType(), TaoADMMRegularizerType, TAOADMM 13636285c0a3SHansol Suh @*/ 13646285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type) 13656285c0a3SHansol Suh { 13666285c0a3SHansol Suh PetscErrorCode ierr; 13676285c0a3SHansol Suh 13686285c0a3SHansol Suh PetscFunctionBegin; 13696285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13706285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao,type,2); 13716285c0a3SHansol Suh ierr = PetscTryMethod(tao,"TaoADMMSetRegularizerType_C",(Tao,TaoADMMRegularizerType),(tao,type));CHKERRQ(ierr); 13726285c0a3SHansol Suh PetscFunctionReturn(0); 13736285c0a3SHansol Suh } 13746285c0a3SHansol Suh 13756285c0a3SHansol Suh /*@ 13766285c0a3SHansol Suh TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM 13776285c0a3SHansol Suh 13786285c0a3SHansol Suh Not Collective 13796285c0a3SHansol Suh 13806285c0a3SHansol Suh Input Parameter: 13816285c0a3SHansol Suh . tao - the Tao context 13826285c0a3SHansol Suh 13836285c0a3SHansol Suh Output Parameter: 13846285c0a3SHansol Suh . type - the type of regularizer 13856285c0a3SHansol Suh 13866285c0a3SHansol Suh Level: intermediate 13876285c0a3SHansol Suh 1388a5a2f7acSBarry Smith .seealso: TaoADMMSetRegularizerType(), TaoADMMRegularizerType, TAOADMM 13896285c0a3SHansol Suh @*/ 13906285c0a3SHansol Suh PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type) 13916285c0a3SHansol Suh { 13926285c0a3SHansol Suh PetscErrorCode ierr; 13936285c0a3SHansol Suh 13946285c0a3SHansol Suh PetscFunctionBegin; 13956285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13966285c0a3SHansol Suh ierr = PetscUseMethod(tao,"TaoADMMGetRegularizerType_C",(Tao,TaoADMMRegularizerType*),(tao,type));CHKERRQ(ierr); 13976285c0a3SHansol Suh PetscFunctionReturn(0); 13986285c0a3SHansol Suh } 13996285c0a3SHansol Suh 14006285c0a3SHansol Suh /*@ 14016285c0a3SHansol Suh TaoADMMSetUpdateType - Set update routine for ADMM routine 14026285c0a3SHansol Suh 14036285c0a3SHansol Suh Not Collective 14046285c0a3SHansol Suh 1405d8d19677SJose E. Roman Input Parameters: 14066285c0a3SHansol Suh + tao - the Tao context 14076285c0a3SHansol Suh - type - spectral parameter update type 14086285c0a3SHansol Suh 14096285c0a3SHansol Suh Level: intermediate 14106285c0a3SHansol Suh 1411a5a2f7acSBarry Smith .seealso: TaoADMMGetUpdateType(), TaoADMMUpdateType, TAOADMM 14126285c0a3SHansol Suh @*/ 14136285c0a3SHansol Suh PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type) 14146285c0a3SHansol Suh { 14156285c0a3SHansol Suh PetscErrorCode ierr; 14166285c0a3SHansol Suh 14176285c0a3SHansol Suh PetscFunctionBegin; 14186285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 14196285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao,type,2); 14206285c0a3SHansol Suh ierr = PetscTryMethod(tao,"TaoADMMSetUpdateType_C",(Tao,TaoADMMUpdateType),(tao,type));CHKERRQ(ierr); 14216285c0a3SHansol Suh PetscFunctionReturn(0); 14226285c0a3SHansol Suh } 14236285c0a3SHansol Suh 14246285c0a3SHansol Suh /*@ 14256285c0a3SHansol Suh TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM 14266285c0a3SHansol Suh 14276285c0a3SHansol Suh Not Collective 14286285c0a3SHansol Suh 14296285c0a3SHansol Suh Input Parameter: 1430f0fc11ceSJed Brown . tao - the Tao context 14316285c0a3SHansol Suh 14326285c0a3SHansol Suh Output Parameter: 14336285c0a3SHansol Suh . type - the type of spectral penalty update routine 14346285c0a3SHansol Suh 14356285c0a3SHansol Suh Level: intermediate 14366285c0a3SHansol Suh 1437a5a2f7acSBarry Smith .seealso: TaoADMMSetUpdateType(), TaoADMMUpdateType, TAOADMM 14386285c0a3SHansol Suh @*/ 14396285c0a3SHansol Suh PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type) 14406285c0a3SHansol Suh { 14416285c0a3SHansol Suh PetscErrorCode ierr; 14426285c0a3SHansol Suh 14436285c0a3SHansol Suh PetscFunctionBegin; 14446285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 14456285c0a3SHansol Suh ierr = PetscUseMethod(tao,"TaoADMMGetUpdateType_C",(Tao,TaoADMMUpdateType*),(tao,type));CHKERRQ(ierr); 14466285c0a3SHansol Suh PetscFunctionReturn(0); 14476285c0a3SHansol Suh } 1448