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 PetscReal Axnorm,Bznorm,ATynorm,temp; 326285c0a3SHansol Suh Vec tempJR,tempL; 336285c0a3SHansol Suh Tao mis; 346285c0a3SHansol Suh 356285c0a3SHansol Suh PetscFunctionBegin; 366285c0a3SHansol Suh mis = am->subsolverX; 376285c0a3SHansol Suh tempJR = am->workJacobianRight; 386285c0a3SHansol Suh tempL = am->workLeft; 396285c0a3SHansol Suh /* ATy */ 409566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(mis, am->y, mis->jacobian_equality, mis->jacobian_equality_pre)); 419566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mis->jacobian_equality,am->y,tempJR)); 429566063dSJacob Faibussowitsch PetscCall(VecNorm(tempJR,NORM_2,&ATynorm)); 436285c0a3SHansol Suh /* dualres = mu * ||AT(Bz-Bzold)||_2 */ 449566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tempJR,-1.,am->Bzold,am->Bz)); 459566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mis->jacobian_equality,tempJR,tempL)); 469566063dSJacob Faibussowitsch PetscCall(VecNorm(tempL,NORM_2,&am->dualres)); 476285c0a3SHansol Suh am->dualres *= am->mu; 486285c0a3SHansol Suh 496285c0a3SHansol Suh /* ||Ax||_2, ||Bz||_2 */ 509566063dSJacob Faibussowitsch PetscCall(VecNorm(am->Ax,NORM_2,&Axnorm)); 519566063dSJacob Faibussowitsch PetscCall(VecNorm(am->Bz,NORM_2,&Bznorm)); 526285c0a3SHansol Suh 536285c0a3SHansol Suh /* Set catol to be catol_admm * max{||Ax||,||Bz||,||c||} * 546285c0a3SHansol Suh * Set gatol to be gatol_admm * ||A^Ty|| * 556285c0a3SHansol Suh * while cnorm is ||r_k||_2, and gnorm is ||d_k||_2 */ 566285c0a3SHansol Suh temp = am->catol_admm * PetscMax(Axnorm, (!am->const_norm) ? Bznorm : PetscMax(Bznorm,am->const_norm)); 579566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintTolerances(tao,temp,PETSC_DEFAULT)); 589566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao, am->gatol_admm*ATynorm, PETSC_DEFAULT,PETSC_DEFAULT)); 596285c0a3SHansol Suh PetscFunctionReturn(0); 606285c0a3SHansol Suh } 616285c0a3SHansol Suh 626285c0a3SHansol Suh /* Penaly Update for Adaptive ADMM. */ 636285c0a3SHansol Suh static PetscErrorCode AdaptiveADMMPenaltyUpdate(Tao tao) 646285c0a3SHansol Suh { 656285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 666285c0a3SHansol Suh PetscReal ydiff_norm, yhatdiff_norm, Axdiff_norm, Bzdiff_norm, Axyhat, Bzy, a_sd, a_mg, a_k, b_sd, b_mg, b_k; 676285c0a3SHansol Suh PetscBool hflag, gflag; 686285c0a3SHansol Suh Vec tempJR,tempJR2; 696285c0a3SHansol Suh 706285c0a3SHansol Suh PetscFunctionBegin; 716285c0a3SHansol Suh tempJR = am->workJacobianRight; 726285c0a3SHansol Suh tempJR2 = am->workJacobianRight2; 736285c0a3SHansol Suh hflag = PETSC_FALSE; 746285c0a3SHansol Suh gflag = PETSC_FALSE; 756285c0a3SHansol Suh a_k = -1; 766285c0a3SHansol Suh b_k = -1; 776285c0a3SHansol Suh 789566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tempJR,-1.,am->Axold,am->Ax)); 799566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tempJR2,-1.,am->yhatold,am->yhat)); 809566063dSJacob Faibussowitsch PetscCall(VecNorm(tempJR,NORM_2,&Axdiff_norm)); 819566063dSJacob Faibussowitsch PetscCall(VecNorm(tempJR2,NORM_2,&yhatdiff_norm)); 829566063dSJacob Faibussowitsch PetscCall(VecDot(tempJR,tempJR2,&Axyhat)); 836285c0a3SHansol Suh 849566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tempJR,-1.,am->Bz0,am->Bz)); 859566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tempJR2,-1.,am->y,am->y0)); 869566063dSJacob Faibussowitsch PetscCall(VecNorm(tempJR,NORM_2,&Bzdiff_norm)); 879566063dSJacob Faibussowitsch PetscCall(VecNorm(tempJR2,NORM_2,&ydiff_norm)); 889566063dSJacob Faibussowitsch PetscCall(VecDot(tempJR,tempJR2,&Bzy)); 896285c0a3SHansol Suh 906285c0a3SHansol Suh if (Axyhat > am->orthval*Axdiff_norm*yhatdiff_norm + am->mueps) { 916285c0a3SHansol Suh hflag = PETSC_TRUE; 926285c0a3SHansol Suh a_sd = PetscSqr(yhatdiff_norm)/Axyhat; /* alphaSD */ 936285c0a3SHansol Suh a_mg = Axyhat/PetscSqr(Axdiff_norm); /* alphaMG */ 946285c0a3SHansol Suh a_k = (a_mg/a_sd) > 0.5 ? a_mg : a_sd - 0.5*a_mg; 956285c0a3SHansol Suh } 966285c0a3SHansol Suh if (Bzy > am->orthval*Bzdiff_norm*ydiff_norm + am->mueps) { 976285c0a3SHansol Suh gflag = PETSC_TRUE; 986285c0a3SHansol Suh b_sd = PetscSqr(ydiff_norm)/Bzy; /* betaSD */ 996285c0a3SHansol Suh b_mg = Bzy/PetscSqr(Bzdiff_norm); /* betaMG */ 1006285c0a3SHansol Suh b_k = (b_mg/b_sd) > 0.5 ? b_mg : b_sd - 0.5*b_mg; 1016285c0a3SHansol Suh } 1026285c0a3SHansol Suh am->muold = am->mu; 1036285c0a3SHansol Suh if (gflag && hflag) { 1046285c0a3SHansol Suh am->mu = PetscSqrtReal(a_k*b_k); 1056285c0a3SHansol Suh } else if (hflag) { 1066285c0a3SHansol Suh am->mu = a_k; 1076285c0a3SHansol Suh } else if (gflag) { 1086285c0a3SHansol Suh am->mu = b_k; 1096285c0a3SHansol Suh } 1106285c0a3SHansol Suh if (am->mu > am->muold) { 1116285c0a3SHansol Suh am->mu = am->muold; 1126285c0a3SHansol Suh } 1136285c0a3SHansol Suh if (am->mu < am->mumin) { 1146285c0a3SHansol Suh am->mu = am->mumin; 1156285c0a3SHansol Suh } 1166285c0a3SHansol Suh PetscFunctionReturn(0); 1176285c0a3SHansol Suh } 1186285c0a3SHansol Suh 1197f5c9be9SBarry Smith static PetscErrorCode TaoADMMSetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType type) 1206285c0a3SHansol Suh { 1216285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1226285c0a3SHansol Suh 1236285c0a3SHansol Suh PetscFunctionBegin; 1246285c0a3SHansol Suh am->regswitch = type; 1256285c0a3SHansol Suh PetscFunctionReturn(0); 1266285c0a3SHansol Suh } 1276285c0a3SHansol Suh 1287f5c9be9SBarry Smith static PetscErrorCode TaoADMMGetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType *type) 1296285c0a3SHansol Suh { 1306285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1316285c0a3SHansol Suh 1326285c0a3SHansol Suh PetscFunctionBegin; 1336285c0a3SHansol Suh *type = am->regswitch; 1346285c0a3SHansol Suh PetscFunctionReturn(0); 1356285c0a3SHansol Suh } 1366285c0a3SHansol Suh 1377f5c9be9SBarry Smith static PetscErrorCode TaoADMMSetUpdateType_ADMM(Tao tao, TaoADMMUpdateType type) 1386285c0a3SHansol Suh { 1396285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1406285c0a3SHansol Suh 1416285c0a3SHansol Suh PetscFunctionBegin; 1426285c0a3SHansol Suh am->update = type; 1436285c0a3SHansol Suh PetscFunctionReturn(0); 1446285c0a3SHansol Suh } 1456285c0a3SHansol Suh 1467f5c9be9SBarry Smith static PetscErrorCode TaoADMMGetUpdateType_ADMM(Tao tao, TaoADMMUpdateType *type) 1476285c0a3SHansol Suh { 1486285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1496285c0a3SHansol Suh 1506285c0a3SHansol Suh PetscFunctionBegin; 1516285c0a3SHansol Suh *type = am->update; 1526285c0a3SHansol Suh PetscFunctionReturn(0); 1536285c0a3SHansol Suh } 1546285c0a3SHansol Suh 1556285c0a3SHansol Suh /* This routine updates Jacobians with new x,z vectors, 1566285c0a3SHansol Suh * and then updates Ax and Bz vectors, then computes updated residual vector*/ 1576285c0a3SHansol Suh static PetscErrorCode ADMMUpdateConstraintResidualVector(Tao tao, Vec x, Vec z, Vec Ax, Vec Bz, Vec residual) 1586285c0a3SHansol Suh { 1596285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1606285c0a3SHansol Suh Tao mis,reg; 1616285c0a3SHansol Suh 1626285c0a3SHansol Suh PetscFunctionBegin; 1636285c0a3SHansol Suh mis = am->subsolverX; 1646285c0a3SHansol Suh reg = am->subsolverZ; 1659566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(mis, x, mis->jacobian_equality, mis->jacobian_equality_pre)); 1669566063dSJacob Faibussowitsch PetscCall(MatMult(mis->jacobian_equality,x,Ax)); 1679566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(reg, z, reg->jacobian_equality, reg->jacobian_equality_pre)); 1689566063dSJacob Faibussowitsch PetscCall(MatMult(reg->jacobian_equality,z,Bz)); 1696285c0a3SHansol Suh 1709566063dSJacob Faibussowitsch PetscCall(VecWAXPY(residual,1.,Bz,Ax)); 1716285c0a3SHansol Suh if (am->constraint != NULL) { 1729566063dSJacob Faibussowitsch PetscCall(VecAXPY(residual,-1.,am->constraint)); 1736285c0a3SHansol Suh } 1746285c0a3SHansol Suh PetscFunctionReturn(0); 1756285c0a3SHansol Suh } 1766285c0a3SHansol Suh 1776285c0a3SHansol Suh /* Updates Augmented Lagrangians to given routines * 1786285c0a3SHansol Suh * For subsolverX, routine needs to be ComputeObjectiveAndGraidnet 1796285c0a3SHansol Suh * Separate Objective and Gradient routines are not supported. */ 1806285c0a3SHansol Suh static PetscErrorCode SubObjGradUpdate(Tao tao, Vec x, PetscReal *f, Vec g, void *ptr) 1816285c0a3SHansol Suh { 1826285c0a3SHansol Suh Tao parent = (Tao)ptr; 1836285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 1846285c0a3SHansol Suh PetscReal temp,temp2; 1856285c0a3SHansol Suh Vec tempJR; 1866285c0a3SHansol Suh 1876285c0a3SHansol Suh PetscFunctionBegin; 1886285c0a3SHansol Suh tempJR = am->workJacobianRight; 1899566063dSJacob Faibussowitsch PetscCall(ADMMUpdateConstraintResidualVector(parent, x, am->subsolverZ->solution, am->Ax, am->Bz, am->residual)); 1909566063dSJacob Faibussowitsch PetscCall((*am->ops->misfitobjgrad)(am->subsolverX,x,f,g,am->misfitobjgradP)); 1916285c0a3SHansol Suh 1926285c0a3SHansol Suh am->last_misfit_val = *f; 1936285c0a3SHansol Suh /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 1949566063dSJacob Faibussowitsch PetscCall(VecTDot(am->residual,am->y,&temp)); 1959566063dSJacob Faibussowitsch PetscCall(VecTDot(am->residual,am->residual,&temp2)); 1966285c0a3SHansol Suh *f += temp + (am->mu/2)*temp2; 1976285c0a3SHansol Suh 1986285c0a3SHansol Suh /* Gradient. Add + mu*AT(Ax+Bz-c) + yTA*/ 1999566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_equality,am->residual,tempJR)); 2009566063dSJacob Faibussowitsch PetscCall(VecAXPY(g,am->mu,tempJR)); 2019566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_equality,am->y,tempJR)); 2029566063dSJacob Faibussowitsch PetscCall(VecAXPY(g,1.,tempJR)); 2036285c0a3SHansol Suh PetscFunctionReturn(0); 2046285c0a3SHansol Suh } 2056285c0a3SHansol Suh 2066285c0a3SHansol Suh /* Updates Augmented Lagrangians to given routines 2076285c0a3SHansol Suh * For subsolverZ, routine needs to be ComputeObjectiveAndGraidnet 2086285c0a3SHansol Suh * Separate Objective and Gradient routines are not supported. */ 2096285c0a3SHansol Suh static PetscErrorCode RegObjGradUpdate(Tao tao, Vec z, PetscReal *f, Vec g, void *ptr) 2106285c0a3SHansol Suh { 2116285c0a3SHansol Suh Tao parent = (Tao)ptr; 2126285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 2136285c0a3SHansol Suh PetscReal temp,temp2; 2146285c0a3SHansol Suh Vec tempJR; 2156285c0a3SHansol Suh 2166285c0a3SHansol Suh PetscFunctionBegin; 2176285c0a3SHansol Suh tempJR = am->workJacobianRight; 2189566063dSJacob Faibussowitsch PetscCall(ADMMUpdateConstraintResidualVector(parent, am->subsolverX->solution, z, am->Ax, am->Bz, am->residual)); 2199566063dSJacob Faibussowitsch PetscCall((*am->ops->regobjgrad)(am->subsolverZ,z,f,g,am->regobjgradP)); 2206285c0a3SHansol Suh am->last_reg_val= *f; 2216285c0a3SHansol Suh /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 2229566063dSJacob Faibussowitsch PetscCall(VecTDot(am->residual,am->y,&temp)); 2239566063dSJacob Faibussowitsch PetscCall(VecTDot(am->residual,am->residual,&temp2)); 2246285c0a3SHansol Suh *f += temp + (am->mu/2)*temp2; 2256285c0a3SHansol Suh 2266285c0a3SHansol Suh /* Gradient. Add + mu*BT(Ax+Bz-c) + yTB*/ 2279566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality,am->residual,tempJR)); 2289566063dSJacob Faibussowitsch PetscCall(VecAXPY(g,am->mu,tempJR)); 2299566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality,am->y,tempJR)); 2309566063dSJacob Faibussowitsch PetscCall(VecAXPY(g,1.,tempJR)); 2316285c0a3SHansol Suh PetscFunctionReturn(0); 2326285c0a3SHansol Suh } 2336285c0a3SHansol Suh 2346285c0a3SHansol Suh /* Computes epsilon padded L1 norm lambda*sum(sqrt(x^2+eps^2)-eps */ 2356285c0a3SHansol Suh static PetscErrorCode ADMML1EpsilonNorm(Tao tao, Vec x, PetscReal eps, PetscReal *norm) 2366285c0a3SHansol Suh { 2376285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 2386285c0a3SHansol Suh PetscInt N; 2396285c0a3SHansol Suh 2406285c0a3SHansol Suh PetscFunctionBegin; 2419566063dSJacob Faibussowitsch PetscCall(VecGetSize(am->workLeft,&N)); 2429566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(am->workLeft,x,x)); 2439566063dSJacob Faibussowitsch PetscCall(VecShift(am->workLeft,am->l1epsilon*am->l1epsilon)); 2449566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(am->workLeft)); 2459566063dSJacob Faibussowitsch PetscCall(VecSum(am->workLeft,norm)); 2466285c0a3SHansol Suh *norm += N*am->l1epsilon; 2476285c0a3SHansol Suh *norm *= am->lambda; 2486285c0a3SHansol Suh PetscFunctionReturn(0); 2496285c0a3SHansol Suh } 2506285c0a3SHansol Suh 2516285c0a3SHansol Suh static PetscErrorCode ADMMInternalHessianUpdate(Mat H, Mat Constraint, PetscBool Identity, void *ptr) 2526285c0a3SHansol Suh { 2536285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)ptr; 2546285c0a3SHansol Suh 2556285c0a3SHansol Suh PetscFunctionBegin; 2566285c0a3SHansol Suh switch (am->update) { 2576285c0a3SHansol Suh case (TAO_ADMM_UPDATE_BASIC): 2586285c0a3SHansol Suh break; 2596285c0a3SHansol Suh case (TAO_ADMM_UPDATE_ADAPTIVE): 2606285c0a3SHansol Suh case (TAO_ADMM_UPDATE_ADAPTIVE_RELAXED): 2616285c0a3SHansol Suh if (H && (am->muold != am->mu)) { 2626285c0a3SHansol Suh if (!Identity) { 2639566063dSJacob Faibussowitsch PetscCall(MatAXPY(H,am->mu-am->muold,Constraint,DIFFERENT_NONZERO_PATTERN)); 2646285c0a3SHansol Suh } else { 2659566063dSJacob Faibussowitsch PetscCall(MatShift(H,am->mu-am->muold)); 2666285c0a3SHansol Suh } 2676285c0a3SHansol Suh } 2686285c0a3SHansol Suh break; 2696285c0a3SHansol Suh } 2706285c0a3SHansol Suh PetscFunctionReturn(0); 2716285c0a3SHansol Suh } 2726285c0a3SHansol Suh 2736285c0a3SHansol Suh /* Updates Hessian - adds second derivative of augmented Lagrangian 2746285c0a3SHansol Suh * H \gets H + \rho*ATA 2756285c0a3SHansol Suh * Here, \rho does not change in TAO_ADMM_UPDATE_BASIC - thus no-op 2766285c0a3SHansol Suh * For ADAPTAIVE,ADAPTIVE_RELAXED, 2776285c0a3SHansol Suh * H \gets H + (\rho-\rhoold)*ATA 2786285c0a3SHansol Suh * Here, we assume that A is linear constraint i.e., doesnt change. 2796285c0a3SHansol Suh * Thus, for both ADAPTIVE, and RELAXED, ATA matrix is pre-set (except for A=I (null case)) see TaoSetUp_ADMM */ 2806285c0a3SHansol Suh static PetscErrorCode SubHessianUpdate(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 2816285c0a3SHansol Suh { 2826285c0a3SHansol Suh Tao parent = (Tao)ptr; 2836285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 2846285c0a3SHansol Suh 2856285c0a3SHansol Suh PetscFunctionBegin; 2866285c0a3SHansol Suh if (am->Hxchange) { 2876285c0a3SHansol Suh /* Case where Hessian gets updated with respect to x vector input. */ 2889566063dSJacob Faibussowitsch PetscCall((*am->ops->misfithess)(am->subsolverX,x,H,Hpre,am->misfithessP)); 2899566063dSJacob Faibussowitsch PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am)); 2906285c0a3SHansol Suh } else if (am->Hxbool) { 2916285c0a3SHansol Suh /* Hessian doesn't get updated. H(x) = c */ 2926285c0a3SHansol Suh /* Update Lagrangian only only per TAO call */ 2939566063dSJacob Faibussowitsch PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am)); 2946285c0a3SHansol Suh am->Hxbool = PETSC_FALSE; 2956285c0a3SHansol Suh } 2966285c0a3SHansol Suh PetscFunctionReturn(0); 2976285c0a3SHansol Suh } 2986285c0a3SHansol Suh 2996285c0a3SHansol Suh /* Same as SubHessianUpdate, except for B matrix instead of A matrix */ 3006285c0a3SHansol Suh static PetscErrorCode RegHessianUpdate(Tao tao, Vec z, Mat H, Mat Hpre, void *ptr) 3016285c0a3SHansol Suh { 3026285c0a3SHansol Suh Tao parent = (Tao)ptr; 3036285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 3046285c0a3SHansol Suh 3056285c0a3SHansol Suh PetscFunctionBegin; 3066285c0a3SHansol Suh 3076285c0a3SHansol Suh if (am->Hzchange) { 3086285c0a3SHansol Suh /* Case where Hessian gets updated with respect to x vector input. */ 3099566063dSJacob Faibussowitsch PetscCall((*am->ops->reghess)(am->subsolverZ,z,H,Hpre,am->reghessP)); 3109566063dSJacob Faibussowitsch PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am)); 3116285c0a3SHansol Suh } else if (am->Hzbool) { 3126285c0a3SHansol Suh /* Hessian doesn't get updated. H(x) = c */ 3136285c0a3SHansol Suh /* Update Lagrangian only only per TAO call */ 3149566063dSJacob Faibussowitsch PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am)); 3156285c0a3SHansol Suh am->Hzbool = PETSC_FALSE; 3166285c0a3SHansol Suh } 3176285c0a3SHansol Suh PetscFunctionReturn(0); 3186285c0a3SHansol Suh } 3196285c0a3SHansol Suh 3206285c0a3SHansol Suh /* Shell Matrix routine for A matrix. 3216285c0a3SHansol Suh * This gets used when user puts NULL for 3226285c0a3SHansol Suh * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...) 3236285c0a3SHansol Suh * Essentially sets A=I*/ 3246285c0a3SHansol Suh static PetscErrorCode JacobianIdentity(Mat mat,Vec in,Vec out) 3256285c0a3SHansol Suh { 3266285c0a3SHansol Suh PetscFunctionBegin; 3279566063dSJacob Faibussowitsch PetscCall(VecCopy(in,out)); 3286285c0a3SHansol Suh PetscFunctionReturn(0); 3296285c0a3SHansol Suh } 3306285c0a3SHansol Suh 3316285c0a3SHansol Suh /* Shell Matrix routine for B matrix. 3326285c0a3SHansol Suh * This gets used when user puts NULL for 3336285c0a3SHansol Suh * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...) 3346285c0a3SHansol Suh * Sets B=-I */ 3356285c0a3SHansol Suh static PetscErrorCode JacobianIdentityB(Mat mat,Vec in,Vec out) 3366285c0a3SHansol Suh { 3376285c0a3SHansol Suh PetscFunctionBegin; 3389566063dSJacob Faibussowitsch PetscCall(VecCopy(in,out)); 3399566063dSJacob Faibussowitsch PetscCall(VecScale(out,-1.)); 3406285c0a3SHansol Suh PetscFunctionReturn(0); 3416285c0a3SHansol Suh } 3426285c0a3SHansol Suh 3436285c0a3SHansol Suh /* Solve f(x) + g(z) s.t. Ax + Bz = c */ 3446285c0a3SHansol Suh static PetscErrorCode TaoSolve_ADMM(Tao tao) 3456285c0a3SHansol Suh { 3466285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 3476285c0a3SHansol Suh PetscInt N; 3486285c0a3SHansol Suh PetscReal reg_func; 3496285c0a3SHansol Suh PetscBool is_reg_shell; 3506285c0a3SHansol Suh Vec tempL; 3516285c0a3SHansol Suh 3526285c0a3SHansol Suh PetscFunctionBegin; 3536285c0a3SHansol Suh if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 3543c859ba3SBarry Smith PetscCheck(am->subsolverX->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetMisfitConstraintJacobian() first"); 3553c859ba3SBarry Smith PetscCheck(am->subsolverZ->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerConstraintJacobian() first"); 3566285c0a3SHansol Suh if (am->constraint != NULL) { 3579566063dSJacob Faibussowitsch PetscCall(VecNorm(am->constraint,NORM_2,&am->const_norm)); 3586285c0a3SHansol Suh } 3596285c0a3SHansol Suh } 3606285c0a3SHansol Suh tempL = am->workLeft; 3619566063dSJacob Faibussowitsch PetscCall(VecGetSize(tempL,&N)); 3626285c0a3SHansol Suh 3636285c0a3SHansol Suh if (am->Hx && am->ops->misfithess) { 3649566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao)); 3656285c0a3SHansol Suh } 3666285c0a3SHansol Suh 3676285c0a3SHansol Suh if (!am->zJI) { 3686285c0a3SHansol Suh /* Currently, B is assumed to be a linear system, i.e., not getting updated*/ 3699566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(am->JB,am->JB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->BTB))); 3706285c0a3SHansol Suh } 3716285c0a3SHansol Suh if (!am->xJI) { 3726285c0a3SHansol Suh /* Currently, A is assumed to be a linear system, i.e., not getting updated*/ 3739566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(am->subsolverX->jacobian_equality,am->subsolverX->jacobian_equality,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->ATA))); 3746285c0a3SHansol Suh } 3756285c0a3SHansol Suh 3766285c0a3SHansol Suh is_reg_shell = PETSC_FALSE; 3776285c0a3SHansol Suh 3789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell)); 3796285c0a3SHansol Suh 3806285c0a3SHansol Suh if (!is_reg_shell) { 3816285c0a3SHansol Suh switch (am->regswitch) { 3826285c0a3SHansol Suh case (TAO_ADMM_REGULARIZER_USER): 3836285c0a3SHansol Suh break; 3846285c0a3SHansol Suh case (TAO_ADMM_REGULARIZER_SOFT_THRESH): 3856285c0a3SHansol Suh /* Soft Threshold. */ 3866285c0a3SHansol Suh break; 3876285c0a3SHansol Suh } 3881baa6e33SBarry Smith if (am->ops->regobjgrad) PetscCall(TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao)); 3896285c0a3SHansol Suh if (am->Hz && am->ops->reghess) { 3909566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao)); 3916285c0a3SHansol Suh } 3926285c0a3SHansol Suh } 3936285c0a3SHansol Suh 3946285c0a3SHansol Suh switch (am->update) { 3956285c0a3SHansol Suh case TAO_ADMM_UPDATE_BASIC: 3966285c0a3SHansol Suh if (am->subsolverX->hessian) { 3976285c0a3SHansol Suh /* In basic case, Hessian does not get updated w.r.t. to spectral penalty 3986285c0a3SHansol Suh * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/ 3996285c0a3SHansol Suh if (!am->xJI) { 4009566063dSJacob Faibussowitsch PetscCall(MatAXPY(am->subsolverX->hessian,am->mu,am->ATA,DIFFERENT_NONZERO_PATTERN)); 4016285c0a3SHansol Suh } else { 4029566063dSJacob Faibussowitsch PetscCall(MatShift(am->subsolverX->hessian,am->mu)); 4036285c0a3SHansol Suh } 4046285c0a3SHansol Suh } 4056285c0a3SHansol Suh if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) { 4066285c0a3SHansol Suh if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) { 4079566063dSJacob Faibussowitsch PetscCall(MatAXPY(am->subsolverZ->hessian,am->mu,am->BTB,DIFFERENT_NONZERO_PATTERN)); 4086285c0a3SHansol Suh } else { 4099566063dSJacob Faibussowitsch PetscCall(MatShift(am->subsolverZ->hessian,am->mu)); 4106285c0a3SHansol Suh } 4116285c0a3SHansol Suh } 4126285c0a3SHansol Suh break; 4136285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE: 4146285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 4156285c0a3SHansol Suh break; 4166285c0a3SHansol Suh } 4176285c0a3SHansol Suh 4189566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation,&cited)); 4196285c0a3SHansol Suh tao->reason = TAO_CONTINUE_ITERATING; 4206285c0a3SHansol Suh 4216285c0a3SHansol Suh while (tao->reason == TAO_CONTINUE_ITERATING) { 422*dbbe0bcdSBarry Smith PetscTryTypeMethod(tao,update, tao->niter, tao->user_update); 4239566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Bz, am->Bzold)); 4246285c0a3SHansol Suh 4256285c0a3SHansol Suh /* x update */ 4269566063dSJacob Faibussowitsch PetscCall(TaoSolve(am->subsolverX)); 4279566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre)); 4289566063dSJacob Faibussowitsch PetscCall(MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution,am->Ax)); 4296285c0a3SHansol Suh 4306285c0a3SHansol Suh am->Hxbool = PETSC_TRUE; 4316285c0a3SHansol Suh 4326285c0a3SHansol Suh /* z update */ 4336285c0a3SHansol Suh switch (am->regswitch) { 4346285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_USER: 4359566063dSJacob Faibussowitsch PetscCall(TaoSolve(am->subsolverZ)); 4366285c0a3SHansol Suh break; 4376285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_SOFT_THRESH: 4386285c0a3SHansol Suh /* L1 assumes A,B jacobians are identity nxn matrix */ 4399566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->workJacobianRight,1/am->mu,am->y,am->Ax)); 4409566063dSJacob Faibussowitsch PetscCall(TaoSoftThreshold(am->workJacobianRight,-am->lambda/am->mu,am->lambda/am->mu,am->subsolverZ->solution)); 4416285c0a3SHansol Suh break; 4426285c0a3SHansol Suh } 4436285c0a3SHansol Suh am->Hzbool = PETSC_TRUE; 4446285c0a3SHansol Suh /* Returns Ax + Bz - c with updated Ax,Bz vectors */ 4459566063dSJacob Faibussowitsch PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual)); 4466285c0a3SHansol Suh /* Dual variable, y += y + mu*(Ax+Bz-c) */ 4479566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->y, am->mu, am->residual, am->yold)); 4486285c0a3SHansol Suh 4496285c0a3SHansol Suh /* stopping tolerance update */ 4509566063dSJacob Faibussowitsch PetscCall(TaoADMMToleranceUpdate(tao)); 4516285c0a3SHansol Suh 4526285c0a3SHansol Suh /* Updating Spectral Penalty */ 4536285c0a3SHansol Suh switch (am->update) { 4546285c0a3SHansol Suh case TAO_ADMM_UPDATE_BASIC: 4556285c0a3SHansol Suh am->muold = am->mu; 4566285c0a3SHansol Suh break; 4576285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE: 4586285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 4596285c0a3SHansol Suh if (tao->niter == 0) { 4609566063dSJacob Faibussowitsch PetscCall(VecCopy(am->y, am->y0)); 4619566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold)); 4621baa6e33SBarry Smith if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint)); 4639566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold)); 4649566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Ax, am->Axold)); 4659566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Bz, am->Bz0)); 4666285c0a3SHansol Suh am->muold = am->mu; 4676285c0a3SHansol Suh } else if (tao->niter % am->T == 1) { 4686285c0a3SHansol Suh /* we have compute Bzold in a previous iteration, and we computed Ax above */ 4699566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold)); 4701baa6e33SBarry Smith if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint)); 4719566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->yhat, -am->mu, am->residual, am->yold)); 4729566063dSJacob Faibussowitsch PetscCall(AdaptiveADMMPenaltyUpdate(tao)); 4739566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Ax, am->Axold)); 4749566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Bz, am->Bz0)); 4759566063dSJacob Faibussowitsch PetscCall(VecCopy(am->yhat, am->yhatold)); 4769566063dSJacob Faibussowitsch PetscCall(VecCopy(am->y, am->y0)); 4776285c0a3SHansol Suh } else { 4786285c0a3SHansol Suh am->muold = am->mu; 4796285c0a3SHansol Suh } 4806285c0a3SHansol Suh break; 4816285c0a3SHansol Suh default: 4826285c0a3SHansol Suh break; 4836285c0a3SHansol Suh } 4846285c0a3SHansol Suh tao->niter++; 4856285c0a3SHansol Suh 4866285c0a3SHansol Suh /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/ 4876285c0a3SHansol Suh switch (am->regswitch) { 4886285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_USER: 4896285c0a3SHansol Suh if (is_reg_shell) { 4909566063dSJacob Faibussowitsch PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func)); 4916285c0a3SHansol Suh } else { 4926285c0a3SHansol Suh (*am->ops->regobjgrad)(am->subsolverZ,am->subsolverX->solution,®_func,tempL,am->regobjgradP); 4936285c0a3SHansol Suh } 4946285c0a3SHansol Suh break; 4956285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_SOFT_THRESH: 4969566063dSJacob Faibussowitsch PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func)); 4976285c0a3SHansol Suh break; 4986285c0a3SHansol Suh } 4999566063dSJacob Faibussowitsch PetscCall(VecCopy(am->y,am->yold)); 5009566063dSJacob Faibussowitsch PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual)); 5019566063dSJacob Faibussowitsch PetscCall(VecNorm(am->residual,NORM_2,&am->resnorm)); 5029566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,am->last_misfit_val + reg_func,am->dualres,am->resnorm,tao->ksp_its)); 5036285c0a3SHansol Suh 5049566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,am->last_misfit_val + reg_func,am->dualres,am->resnorm,1.0)); 505*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 5066285c0a3SHansol Suh } 5076285c0a3SHansol Suh /* Update vectors */ 5089566063dSJacob Faibussowitsch PetscCall(VecCopy(am->subsolverX->solution,tao->solution)); 5099566063dSJacob Faibussowitsch PetscCall(VecCopy(am->subsolverX->gradient,tao->gradient)); 5109566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", NULL)); 5119566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", NULL)); 5129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL)); 5139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL)); 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL)); 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL)); 5166285c0a3SHansol Suh PetscFunctionReturn(0); 5176285c0a3SHansol Suh } 5186285c0a3SHansol Suh 519*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_ADMM(Tao tao,PetscOptionItems *PetscOptionsObject) 5206285c0a3SHansol Suh { 5216285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5226285c0a3SHansol Suh 5236285c0a3SHansol Suh PetscFunctionBegin; 524d0609cedSBarry Smith PetscOptionsHeadBegin(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. "); 5259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient","regularizer constant","",am->lambda,&am->lambda,NULL)); 5269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty","Constant for Augmented Lagrangian term.","",am->mu,&am->mu,NULL)); 5279566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter","x relaxation parameter for Z update.","",am->gamma,&am->gamma,NULL)); 5289566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor","ADMM dynamic tolerance update factor.","",am->tol,&am->tol,NULL)); 5299566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor","ADMM spectral penalty update curvature safeguard value.","",am->orthval,&am->orthval,NULL)); 5309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty","Set ADMM minimum spectral penalty.","",am->mumin,&am->mumin,NULL)); 531d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-tao_admm_dual_update","Lagrangian dual update policy","TaoADMMUpdateType",TaoADMMUpdateTypes,(PetscEnum)am->update,(PetscEnum*)&am->update,NULL)); 532d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-tao_admm_regularizer_type","ADMM regularizer update rule","TaoADMMRegularizerType",TaoADMMRegularizerTypes,(PetscEnum)am->regswitch,(PetscEnum*)&am->regswitch,NULL)); 533d0609cedSBarry Smith PetscOptionsHeadEnd(); 5349566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(am->subsolverX)); 5356285c0a3SHansol Suh if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 5369566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(am->subsolverZ)); 5376285c0a3SHansol Suh } 5386285c0a3SHansol Suh PetscFunctionReturn(0); 5396285c0a3SHansol Suh } 5406285c0a3SHansol Suh 5416285c0a3SHansol Suh static PetscErrorCode TaoView_ADMM(Tao tao,PetscViewer viewer) 5426285c0a3SHansol Suh { 5436285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5446285c0a3SHansol Suh 5456285c0a3SHansol Suh PetscFunctionBegin; 5469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 5479566063dSJacob Faibussowitsch PetscCall(TaoView(am->subsolverX,viewer)); 5489566063dSJacob Faibussowitsch PetscCall(TaoView(am->subsolverZ,viewer)); 5499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 5506285c0a3SHansol Suh PetscFunctionReturn(0); 5516285c0a3SHansol Suh } 5526285c0a3SHansol Suh 5536285c0a3SHansol Suh static PetscErrorCode TaoSetUp_ADMM(Tao tao) 5546285c0a3SHansol Suh { 5556285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5566285c0a3SHansol Suh PetscInt n,N,M; 5576285c0a3SHansol Suh 5586285c0a3SHansol Suh PetscFunctionBegin; 5599566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution,&n)); 5609566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution,&N)); 5616285c0a3SHansol Suh /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */ 5626285c0a3SHansol Suh if (!am->JB) { 5636285c0a3SHansol Suh am->zJI = PETSC_TRUE; 5649566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JB)); 5659566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(am->JB,MATOP_MULT,(void (*)(void))JacobianIdentityB)); 5669566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(am->JB,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentityB)); 5676285c0a3SHansol Suh am->JBpre = am->JB; 5686285c0a3SHansol Suh } 5696285c0a3SHansol Suh if (!am->JA) { 5706285c0a3SHansol Suh am->xJI = PETSC_TRUE; 5719566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JA)); 5729566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(am->JA,MATOP_MULT,(void (*)(void))JacobianIdentity)); 5739566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(am->JA,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentity)); 5746285c0a3SHansol Suh am->JApre = am->JA; 5756285c0a3SHansol Suh } 5769566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(am->JA,NULL,&am->Ax)); 5776285c0a3SHansol Suh if (!tao->gradient) { 5789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 5796285c0a3SHansol Suh } 5809566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(am->subsolverX, tao->solution)); 5816285c0a3SHansol Suh if (!am->z) { 5829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&am->z)); 5839566063dSJacob Faibussowitsch PetscCall(VecSet(am->z,0.0)); 5846285c0a3SHansol Suh } 5859566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(am->subsolverZ, am->z)); 5866285c0a3SHansol Suh if (!am->workLeft) { 5879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&am->workLeft)); 5886285c0a3SHansol Suh } 5896285c0a3SHansol Suh if (!am->Axold) { 5909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->Axold)); 5916285c0a3SHansol Suh } 5926285c0a3SHansol Suh if (!am->workJacobianRight) { 5939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight)); 5946285c0a3SHansol Suh } 5956285c0a3SHansol Suh if (!am->workJacobianRight2) { 5969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight2)); 5976285c0a3SHansol Suh } 5986285c0a3SHansol Suh if (!am->Bz) { 5999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->Bz)); 6006285c0a3SHansol Suh } 6016285c0a3SHansol Suh if (!am->Bzold) { 6029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->Bzold)); 6036285c0a3SHansol Suh } 6046285c0a3SHansol Suh if (!am->Bz0) { 6059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->Bz0)); 6066285c0a3SHansol Suh } 6076285c0a3SHansol Suh if (!am->y) { 6089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->y)); 6099566063dSJacob Faibussowitsch PetscCall(VecSet(am->y,0.0)); 6106285c0a3SHansol Suh } 6116285c0a3SHansol Suh if (!am->yold) { 6129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->yold)); 6139566063dSJacob Faibussowitsch PetscCall(VecSet(am->yold,0.0)); 6146285c0a3SHansol Suh } 6156285c0a3SHansol Suh if (!am->y0) { 6169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->y0)); 6179566063dSJacob Faibussowitsch PetscCall(VecSet(am->y0,0.0)); 6186285c0a3SHansol Suh } 6196285c0a3SHansol Suh if (!am->yhat) { 6209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->yhat)); 6219566063dSJacob Faibussowitsch PetscCall(VecSet(am->yhat,0.0)); 6226285c0a3SHansol Suh } 6236285c0a3SHansol Suh if (!am->yhatold) { 6249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->yhatold)); 6259566063dSJacob Faibussowitsch PetscCall(VecSet(am->yhatold,0.0)); 6266285c0a3SHansol Suh } 6276285c0a3SHansol Suh if (!am->residual) { 6289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax,&am->residual)); 6299566063dSJacob Faibussowitsch PetscCall(VecSet(am->residual,0.0)); 6306285c0a3SHansol Suh } 6316285c0a3SHansol Suh if (!am->constraint) { 6326285c0a3SHansol Suh am->constraint = NULL; 6336285c0a3SHansol Suh } else { 6349566063dSJacob Faibussowitsch PetscCall(VecGetSize(am->constraint,&M)); 6353c859ba3SBarry Smith PetscCheck(M == N,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Solution vector and constraint vector must be of same size!"); 6366285c0a3SHansol Suh } 6376285c0a3SHansol Suh 6386285c0a3SHansol Suh /* Save changed tao tolerance for adaptive tolerance */ 6396285c0a3SHansol Suh if (tao->gatol_changed) { 6406285c0a3SHansol Suh am->gatol_admm = tao->gatol; 6416285c0a3SHansol Suh } 6426285c0a3SHansol Suh if (tao->catol_changed) { 6436285c0a3SHansol Suh am->catol_admm = tao->catol; 6446285c0a3SHansol Suh } 6456285c0a3SHansol Suh 6466285c0a3SHansol Suh /*Update spectral and dual elements to X subsolver */ 6479566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao)); 6489566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX,am->JA,am->JApre, am->ops->misfitjac, am->misfitjacobianP)); 6499566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ,am->JB,am->JBpre, am->ops->regjac, am->regjacobianP)); 6506285c0a3SHansol Suh PetscFunctionReturn(0); 6516285c0a3SHansol Suh } 6526285c0a3SHansol Suh 6536285c0a3SHansol Suh static PetscErrorCode TaoDestroy_ADMM(Tao tao) 6546285c0a3SHansol Suh { 6556285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 6566285c0a3SHansol Suh 6576285c0a3SHansol Suh PetscFunctionBegin; 6589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->z)); 6599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Ax)); 6609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Axold)); 6619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Bz)); 6629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Bzold)); 6639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Bz0)); 6649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->residual)); 6659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->y)); 6669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->yold)); 6679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->y0)); 6689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->yhat)); 6699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->yhatold)); 6709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->workLeft)); 6719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->workJacobianRight)); 6729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->workJacobianRight2)); 6736285c0a3SHansol Suh 6749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JA)); 6759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JB)); 6766285c0a3SHansol Suh if (!am->xJI) { 6779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JApre)); 6786285c0a3SHansol Suh } 6796285c0a3SHansol Suh if (!am->zJI) { 6809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JBpre)); 6816285c0a3SHansol Suh } 6826285c0a3SHansol Suh if (am->Hx) { 6839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hx)); 6849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hxpre)); 6856285c0a3SHansol Suh } 6866285c0a3SHansol Suh if (am->Hz) { 6879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hz)); 6889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hzpre)); 6896285c0a3SHansol Suh } 6909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->ATA)); 6919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->BTB)); 6929566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&am->subsolverX)); 6939566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&am->subsolverZ)); 6946285c0a3SHansol Suh am->parent = NULL; 6952e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL)); 6962e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL)); 6972e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL)); 6982e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL)); 6999566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 7006285c0a3SHansol Suh PetscFunctionReturn(0); 7016285c0a3SHansol Suh } 7026285c0a3SHansol Suh 7036285c0a3SHansol Suh /*MC 7046285c0a3SHansol Suh 7056285c0a3SHansol Suh TAOADMM - Alternating direction method of multipliers method fo solving linear problems with 7066285c0a3SHansol Suh constraints. in a min_x f(x) + g(z) s.t. Ax+Bz=c. 707a5b23f4aSJose E. Roman This algorithm employs two sub Tao solvers, of which type can be specified 7086285c0a3SHansol Suh by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers. 7096285c0a3SHansol Suh Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via 7106285c0a3SHansol Suh TaoADMMSet{Misfit,Regularizer}HessianChangeStatus. 7116285c0a3SHansol Suh Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type. 7126285c0a3SHansol Suh There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update, 713a5b23f4aSJose E. Roman currently there are basic option and adaptive option. 7146285c0a3SHansol Suh Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian. 7156285c0a3SHansol Suh c can be set with TaoADMMSetConstraintVectorRHS. 7166285c0a3SHansol Suh The user can also provide regularizer weight for second subsolver. 7176285c0a3SHansol Suh 7186285c0a3SHansol Suh References: 719606c0280SSatish Balay . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom 7206285c0a3SHansol Suh "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation" 7216285c0a3SHansol Suh The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017. 7226285c0a3SHansol Suh 7236285c0a3SHansol Suh Options Database Keys: 7246285c0a3SHansol Suh + -tao_admm_regularizer_coefficient - regularizer constant (default 1.e-6) 7256285c0a3SHansol Suh . -tao_admm_spectral_penalty - Constant for Augmented Lagrangian term (default 1.) 7266285c0a3SHansol Suh . -tao_admm_relaxation_parameter - relaxation parameter for Z update (default 1.) 7276285c0a3SHansol Suh . -tao_admm_tolerance_update_factor - ADMM dynamic tolerance update factor (default 1.e-12) 7286285c0a3SHansol Suh . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2) 7296285c0a3SHansol Suh . -tao_admm_minimum_spectral_penalty - Set ADMM minimum spectral penalty (default 0) 7306285c0a3SHansol Suh . -tao_admm_dual_update - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic") 7316285c0a3SHansol Suh - -tao_admm_regularizer_type - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold") 7326285c0a3SHansol Suh 7336285c0a3SHansol Suh Level: beginner 7346285c0a3SHansol Suh 735db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitHessianChangeStatus()`, `TaoADMMSetRegHessianChangeStatus()`, `TaoADMMGetSpectralPenalty()`, 736db781477SPatrick Sanan `TaoADMMGetMisfitSubsolver()`, `TaoADMMGetRegularizationSubsolver()`, `TaoADMMSetConstraintVectorRHS()`, 737db781477SPatrick Sanan `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetRegularizerCoefficient()`, 738db781477SPatrick Sanan `TaoADMMSetRegularizerConstraintJacobian()`, `TaoADMMSetMisfitConstraintJacobian()`, 739db781477SPatrick Sanan `TaoADMMSetMisfitObjectiveAndGradientRoutine()`, `TaoADMMSetMisfitHessianRoutine()`, 740db781477SPatrick Sanan `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, `TaoADMMSetRegularizerHessianRoutine()`, 741db781477SPatrick Sanan `TaoGetADMMParentTao()`, `TaoADMMGetDualVector()`, `TaoADMMSetRegularizerType()`, 742db781477SPatrick Sanan `TaoADMMGetRegularizerType()`, `TaoADMMSetUpdateType()`, `TaoADMMGetUpdateType()` 7436285c0a3SHansol Suh M*/ 7446285c0a3SHansol Suh 7456285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao) 7466285c0a3SHansol Suh { 7476285c0a3SHansol Suh TAO_ADMM *am; 7486285c0a3SHansol Suh 7496285c0a3SHansol Suh PetscFunctionBegin; 7509566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&am)); 7516285c0a3SHansol Suh 7526285c0a3SHansol Suh tao->ops->destroy = TaoDestroy_ADMM; 7536285c0a3SHansol Suh tao->ops->setup = TaoSetUp_ADMM; 7546285c0a3SHansol Suh tao->ops->setfromoptions = TaoSetFromOptions_ADMM; 7556285c0a3SHansol Suh tao->ops->view = TaoView_ADMM; 7566285c0a3SHansol Suh tao->ops->solve = TaoSolve_ADMM; 7576285c0a3SHansol Suh 7586285c0a3SHansol Suh tao->data = (void*)am; 7596285c0a3SHansol Suh am->l1epsilon = 1e-6; 7606285c0a3SHansol Suh am->lambda = 1e-4; 7616285c0a3SHansol Suh am->mu = 1.; 7626285c0a3SHansol Suh am->muold = 0.; 7636285c0a3SHansol Suh am->mueps = PETSC_MACHINE_EPSILON; 7646285c0a3SHansol Suh am->mumin = 0.; 7656285c0a3SHansol Suh am->orthval = 0.2; 7666285c0a3SHansol Suh am->T = 2; 7676285c0a3SHansol Suh am->parent = tao; 7686285c0a3SHansol Suh am->update = TAO_ADMM_UPDATE_BASIC; 7696285c0a3SHansol Suh am->regswitch = TAO_ADMM_REGULARIZER_SOFT_THRESH; 7706285c0a3SHansol Suh am->tol = PETSC_SMALL; 7716285c0a3SHansol Suh am->const_norm = 0; 7726285c0a3SHansol Suh am->resnorm = 0; 7736285c0a3SHansol Suh am->dualres = 0; 77483c8fe1dSLisandro Dalcin am->ops->regobjgrad = NULL; 77583c8fe1dSLisandro Dalcin am->ops->reghess = NULL; 7766285c0a3SHansol Suh am->gamma = 1; 77783c8fe1dSLisandro Dalcin am->regobjgradP = NULL; 77883c8fe1dSLisandro Dalcin am->reghessP = NULL; 7796285c0a3SHansol Suh am->gatol_admm = 1e-8; 7806285c0a3SHansol Suh am->catol_admm = 0; 7816285c0a3SHansol Suh am->Hxchange = PETSC_TRUE; 7826285c0a3SHansol Suh am->Hzchange = PETSC_TRUE; 7836285c0a3SHansol Suh am->Hzbool = PETSC_TRUE; 7846285c0a3SHansol Suh am->Hxbool = PETSC_TRUE; 7856285c0a3SHansol Suh 7869566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverX)); 7879566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(am->subsolverX,"misfit_")); 7889566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX,(PetscObject)tao,1)); 7899566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverZ)); 7909566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(am->subsolverZ,"reg_")); 7919566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ,(PetscObject)tao,1)); 7926285c0a3SHansol Suh 7939566063dSJacob Faibussowitsch PetscCall(TaoSetType(am->subsolverX,TAONLS)); 7949566063dSJacob Faibussowitsch PetscCall(TaoSetType(am->subsolverZ,TAONLS)); 7959566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", (PetscObject) tao)); 7969566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", (PetscObject) tao)); 7979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",TaoADMMSetRegularizerType_ADMM)); 7989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",TaoADMMGetRegularizerType_ADMM)); 7999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",TaoADMMSetUpdateType_ADMM)); 8009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",TaoADMMGetUpdateType_ADMM)); 8016285c0a3SHansol Suh PetscFunctionReturn(0); 8026285c0a3SHansol Suh } 8036285c0a3SHansol Suh 8046285c0a3SHansol Suh /*@ 8056285c0a3SHansol Suh TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines whether Hessian matrix of misfit subsolver changes with respect to input vector. 8066285c0a3SHansol Suh 8076285c0a3SHansol Suh Collective on Tao 8086285c0a3SHansol Suh 8096285c0a3SHansol Suh Input Parameters: 8107f5c9be9SBarry Smith + tao - the Tao solver context. 8117f5c9be9SBarry Smith - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 8126285c0a3SHansol Suh 8136285c0a3SHansol Suh Level: advanced 814a5a2f7acSBarry Smith 815db781477SPatrick Sanan .seealso: `TAOADMM` 816a5a2f7acSBarry Smith 8176285c0a3SHansol Suh @*/ 8186285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b) 8196285c0a3SHansol Suh { 8206285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8216285c0a3SHansol Suh 8226285c0a3SHansol Suh PetscFunctionBegin; 8236285c0a3SHansol Suh am->Hxchange = b; 8246285c0a3SHansol Suh PetscFunctionReturn(0); 8256285c0a3SHansol Suh } 8266285c0a3SHansol Suh 8276285c0a3SHansol Suh /*@ 8286285c0a3SHansol Suh TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector. 8297f5c9be9SBarry Smith 8306285c0a3SHansol Suh Collective on Tao 8316285c0a3SHansol Suh 8326285c0a3SHansol Suh Input Parameters: 8336285c0a3SHansol Suh + tao - the Tao solver context 8347f5c9be9SBarry Smith - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 8356285c0a3SHansol Suh 8366285c0a3SHansol Suh Level: advanced 837a5a2f7acSBarry Smith 838db781477SPatrick Sanan .seealso: `TAOADMM` 839a5a2f7acSBarry Smith 8406285c0a3SHansol Suh @*/ 8416285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b) 8426285c0a3SHansol Suh { 8436285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8446285c0a3SHansol Suh 8456285c0a3SHansol Suh PetscFunctionBegin; 8466285c0a3SHansol Suh am->Hzchange = b; 8476285c0a3SHansol Suh PetscFunctionReturn(0); 8486285c0a3SHansol Suh } 8496285c0a3SHansol Suh 8506285c0a3SHansol Suh /*@ 8516285c0a3SHansol Suh TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value 8526285c0a3SHansol Suh 8536285c0a3SHansol Suh Collective on Tao 8546285c0a3SHansol Suh 8556285c0a3SHansol Suh Input Parameters: 8566285c0a3SHansol Suh + tao - the Tao solver context 8577f5c9be9SBarry Smith - mu - spectral penalty 8586285c0a3SHansol Suh 8596285c0a3SHansol Suh Level: advanced 8606285c0a3SHansol Suh 861db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TAOADMM` 8626285c0a3SHansol Suh @*/ 8636285c0a3SHansol Suh PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu) 8646285c0a3SHansol Suh { 8656285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8666285c0a3SHansol Suh 8676285c0a3SHansol Suh PetscFunctionBegin; 8686285c0a3SHansol Suh am->mu = mu; 8696285c0a3SHansol Suh PetscFunctionReturn(0); 8706285c0a3SHansol Suh } 8716285c0a3SHansol Suh 8726285c0a3SHansol Suh /*@ 8736285c0a3SHansol Suh TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value 8746285c0a3SHansol Suh 8756285c0a3SHansol Suh Collective on Tao 8766285c0a3SHansol Suh 8777f5c9be9SBarry Smith Input Parameter: 8787f5c9be9SBarry Smith . tao - the Tao solver context 8797f5c9be9SBarry Smith 8807f5c9be9SBarry Smith Output Parameter: 8817f5c9be9SBarry Smith . mu - spectral penalty 8826285c0a3SHansol Suh 8836285c0a3SHansol Suh Level: advanced 8846285c0a3SHansol Suh 885db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetSpectralPenalty()`, `TAOADMM` 8866285c0a3SHansol Suh @*/ 8876285c0a3SHansol Suh PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu) 8886285c0a3SHansol Suh { 8896285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8906285c0a3SHansol Suh 8916285c0a3SHansol Suh PetscFunctionBegin; 8926285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 8936285c0a3SHansol Suh PetscValidRealPointer(mu,2); 8946285c0a3SHansol Suh *mu = am->mu; 8956285c0a3SHansol Suh PetscFunctionReturn(0); 8966285c0a3SHansol Suh } 8976285c0a3SHansol Suh 8986285c0a3SHansol Suh /*@ 8996285c0a3SHansol Suh TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM 9006285c0a3SHansol Suh 9016285c0a3SHansol Suh Collective on Tao 9026285c0a3SHansol Suh 9037f5c9be9SBarry Smith Input Parameter: 9047f5c9be9SBarry Smith . tao - the Tao solver context 9057f5c9be9SBarry Smith 9067f5c9be9SBarry Smith Output Parameter: 9077f5c9be9SBarry Smith . misfit - the Tao subsolver context 9086285c0a3SHansol Suh 9096285c0a3SHansol Suh Level: advanced 910a5a2f7acSBarry Smith 911db781477SPatrick Sanan .seealso: `TAOADMM` 912a5a2f7acSBarry Smith 9136285c0a3SHansol Suh @*/ 9146285c0a3SHansol Suh PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit) 9156285c0a3SHansol Suh { 9166285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9176285c0a3SHansol Suh 9186285c0a3SHansol Suh PetscFunctionBegin; 9196285c0a3SHansol Suh *misfit = am->subsolverX; 9206285c0a3SHansol Suh PetscFunctionReturn(0); 9216285c0a3SHansol Suh } 9226285c0a3SHansol Suh 9236285c0a3SHansol Suh /*@ 9246285c0a3SHansol Suh TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM 9256285c0a3SHansol Suh 9266285c0a3SHansol Suh Collective on Tao 9276285c0a3SHansol Suh 9287f5c9be9SBarry Smith Input Parameter: 9297f5c9be9SBarry Smith . tao - the Tao solver context 9307f5c9be9SBarry Smith 9317f5c9be9SBarry Smith Output Parameter: 9327f5c9be9SBarry Smith . reg - the Tao subsolver context 9336285c0a3SHansol Suh 9346285c0a3SHansol Suh Level: advanced 935a5a2f7acSBarry Smith 936db781477SPatrick Sanan .seealso: `TAOADMM` 937a5a2f7acSBarry Smith 9386285c0a3SHansol Suh @*/ 9396285c0a3SHansol Suh PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg) 9406285c0a3SHansol Suh { 9416285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9426285c0a3SHansol Suh 9436285c0a3SHansol Suh PetscFunctionBegin; 9446285c0a3SHansol Suh *reg = am->subsolverZ; 9456285c0a3SHansol Suh PetscFunctionReturn(0); 9466285c0a3SHansol Suh } 9476285c0a3SHansol Suh 9486285c0a3SHansol Suh /*@ 9496285c0a3SHansol Suh TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM 9506285c0a3SHansol Suh 9516285c0a3SHansol Suh Collective on Tao 9526285c0a3SHansol Suh 9536285c0a3SHansol Suh Input Parameters: 9546285c0a3SHansol Suh + tao - the Tao solver context 9556285c0a3SHansol Suh - c - RHS vector 9566285c0a3SHansol Suh 9576285c0a3SHansol Suh Level: advanced 958a5a2f7acSBarry Smith 959db781477SPatrick Sanan .seealso: `TAOADMM` 960a5a2f7acSBarry Smith 9616285c0a3SHansol Suh @*/ 9626285c0a3SHansol Suh PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c) 9636285c0a3SHansol Suh { 9646285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9656285c0a3SHansol Suh 9666285c0a3SHansol Suh PetscFunctionBegin; 9676285c0a3SHansol Suh am->constraint = c; 9686285c0a3SHansol Suh PetscFunctionReturn(0); 9696285c0a3SHansol Suh } 9706285c0a3SHansol Suh 9716285c0a3SHansol Suh /*@ 9727f5c9be9SBarry Smith TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty 9736285c0a3SHansol Suh 9746285c0a3SHansol Suh Collective on Tao 9756285c0a3SHansol Suh 9766285c0a3SHansol Suh Input Parameters: 9776285c0a3SHansol Suh + tao - the Tao solver context 9786285c0a3SHansol Suh - mu - minimum spectral penalty value 9796285c0a3SHansol Suh 9806285c0a3SHansol Suh Level: advanced 9816285c0a3SHansol Suh 982db781477SPatrick Sanan .seealso: `TaoADMMGetSpectralPenalty()`, `TAOADMM` 9836285c0a3SHansol Suh @*/ 9846285c0a3SHansol Suh PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu) 9856285c0a3SHansol Suh { 9866285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9876285c0a3SHansol Suh 9886285c0a3SHansol Suh PetscFunctionBegin; 9896285c0a3SHansol Suh am->mumin= mu; 9906285c0a3SHansol Suh PetscFunctionReturn(0); 9916285c0a3SHansol Suh } 9926285c0a3SHansol Suh 9936285c0a3SHansol Suh /*@ 9946285c0a3SHansol Suh TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case 9956285c0a3SHansol Suh 9966285c0a3SHansol Suh Collective on Tao 9976285c0a3SHansol Suh 9986285c0a3SHansol Suh Input Parameters: 9996285c0a3SHansol Suh + tao - the Tao solver context 10006285c0a3SHansol Suh - lambda - L1-norm regularizer coefficient 10016285c0a3SHansol Suh 10026285c0a3SHansol Suh Level: advanced 10037f5c9be9SBarry Smith 1004db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM` 10057f5c9be9SBarry Smith 10066285c0a3SHansol Suh @*/ 10076285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda) 10086285c0a3SHansol Suh { 10096285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10106285c0a3SHansol Suh 10116285c0a3SHansol Suh PetscFunctionBegin; 10126285c0a3SHansol Suh am->lambda = lambda; 10136285c0a3SHansol Suh PetscFunctionReturn(0); 10146285c0a3SHansol Suh } 10156285c0a3SHansol Suh 10166285c0a3SHansol Suh /*@C 10177f5c9be9SBarry Smith TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable. 10186285c0a3SHansol Suh 10196285c0a3SHansol Suh Collective on Tao 10206285c0a3SHansol Suh 1021a5a2f7acSBarry Smith Input Parameters: 1022a5a2f7acSBarry Smith + tao - the Tao solver context 10236285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix 10246285c0a3SHansol Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 10256285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function 10266285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 10276285c0a3SHansol Suh 10286285c0a3SHansol Suh Level: advanced 10297f5c9be9SBarry Smith 1030db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM` 10317f5c9be9SBarry Smith 10326285c0a3SHansol Suh @*/ 10336285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 10346285c0a3SHansol Suh { 10356285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10366285c0a3SHansol Suh 10376285c0a3SHansol Suh PetscFunctionBegin; 10386285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 10396285c0a3SHansol Suh if (J) { 10406285c0a3SHansol Suh PetscValidHeaderSpecific(J,MAT_CLASSID,2); 10416285c0a3SHansol Suh PetscCheckSameComm(tao,1,J,2); 10426285c0a3SHansol Suh } 10436285c0a3SHansol Suh if (Jpre) { 10446285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 10456285c0a3SHansol Suh PetscCheckSameComm(tao,1,Jpre,3); 10466285c0a3SHansol Suh } 10476285c0a3SHansol Suh if (ctx) am->misfitjacobianP = ctx; 10486285c0a3SHansol Suh if (func) am->ops->misfitjac = func; 10496285c0a3SHansol Suh 10506285c0a3SHansol Suh if (J) { 10519566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J)); 10529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JA)); 10536285c0a3SHansol Suh am->JA = J; 10546285c0a3SHansol Suh } 10556285c0a3SHansol Suh if (Jpre) { 10569566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre)); 10579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JApre)); 10586285c0a3SHansol Suh am->JApre = Jpre; 10596285c0a3SHansol Suh } 10606285c0a3SHansol Suh PetscFunctionReturn(0); 10616285c0a3SHansol Suh } 10626285c0a3SHansol Suh 10636285c0a3SHansol Suh /*@C 10646285c0a3SHansol Suh TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable. 10656285c0a3SHansol Suh 10666285c0a3SHansol Suh Collective on Tao 10676285c0a3SHansol Suh 10686285c0a3SHansol Suh Input Parameters: 10696285c0a3SHansol Suh + tao - the Tao solver context 10706285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix 10716285c0a3SHansol Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 10726285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function 10736285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 10746285c0a3SHansol Suh 10756285c0a3SHansol Suh Level: advanced 10767f5c9be9SBarry Smith 1077db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetMisfitConstraintJacobian()`, `TAOADMM` 10787f5c9be9SBarry Smith 10796285c0a3SHansol Suh @*/ 10806285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 10816285c0a3SHansol Suh { 10826285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10836285c0a3SHansol Suh 10846285c0a3SHansol Suh PetscFunctionBegin; 10856285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 10866285c0a3SHansol Suh if (J) { 10876285c0a3SHansol Suh PetscValidHeaderSpecific(J,MAT_CLASSID,2); 10886285c0a3SHansol Suh PetscCheckSameComm(tao,1,J,2); 10896285c0a3SHansol Suh } 10906285c0a3SHansol Suh if (Jpre) { 10916285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 10926285c0a3SHansol Suh PetscCheckSameComm(tao,1,Jpre,3); 10936285c0a3SHansol Suh } 10946285c0a3SHansol Suh if (ctx) am->regjacobianP = ctx; 10956285c0a3SHansol Suh if (func) am->ops->regjac = func; 10966285c0a3SHansol Suh 10976285c0a3SHansol Suh if (J) { 10989566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J)); 10999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JB)); 11006285c0a3SHansol Suh am->JB = J; 11016285c0a3SHansol Suh } 11026285c0a3SHansol Suh if (Jpre) { 11039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre)); 11049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JBpre)); 11056285c0a3SHansol Suh am->JBpre = Jpre; 11066285c0a3SHansol Suh } 11076285c0a3SHansol Suh PetscFunctionReturn(0); 11086285c0a3SHansol Suh } 11096285c0a3SHansol Suh 11106285c0a3SHansol Suh /*@C 11117f5c9be9SBarry Smith TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function 11127f5c9be9SBarry Smith 11137f5c9be9SBarry Smith Collective on tao 11146285c0a3SHansol Suh 11156285c0a3SHansol Suh Input Parameters: 11166285c0a3SHansol Suh + tao - the Tao context 11176285c0a3SHansol Suh . func - function pointer for the misfit value and gradient evaluation 11186285c0a3SHansol Suh - ctx - user context for the misfit 11196285c0a3SHansol Suh 11206285c0a3SHansol Suh Level: advanced 11217f5c9be9SBarry Smith 1122db781477SPatrick Sanan .seealso: `TAOADMM` 11237f5c9be9SBarry Smith 11246285c0a3SHansol Suh @*/ 11256285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 11266285c0a3SHansol Suh { 11276285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11286285c0a3SHansol Suh 11296285c0a3SHansol Suh PetscFunctionBegin; 11306285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11316285c0a3SHansol Suh am->misfitobjgradP = ctx; 11326285c0a3SHansol Suh am->ops->misfitobjgrad = func; 11336285c0a3SHansol Suh PetscFunctionReturn(0); 11346285c0a3SHansol Suh } 11356285c0a3SHansol Suh 11366285c0a3SHansol Suh /*@C 11376285c0a3SHansol Suh TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back 11386285c0a3SHansol Suh function into the algorithm, to be used for subsolverX. 11396285c0a3SHansol Suh 11407f5c9be9SBarry Smith Collective on tao 11417f5c9be9SBarry Smith 11426285c0a3SHansol Suh Input Parameters: 11436285c0a3SHansol Suh + tao - the Tao context 11446285c0a3SHansol Suh . H - user-created matrix for the Hessian of the misfit term 11456285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term 11466285c0a3SHansol Suh . func - function pointer for the misfit Hessian evaluation 11476285c0a3SHansol Suh - ctx - user context for the misfit Hessian 11486285c0a3SHansol Suh 11496285c0a3SHansol Suh Level: advanced 1150a5a2f7acSBarry Smith 1151db781477SPatrick Sanan .seealso: `TAOADMM` 1152a5a2f7acSBarry Smith 11536285c0a3SHansol Suh @*/ 11546285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 11556285c0a3SHansol Suh { 11566285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11576285c0a3SHansol Suh 11586285c0a3SHansol Suh PetscFunctionBegin; 11596285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11606285c0a3SHansol Suh if (H) { 11616285c0a3SHansol Suh PetscValidHeaderSpecific(H,MAT_CLASSID,2); 11626285c0a3SHansol Suh PetscCheckSameComm(tao,1,H,2); 11636285c0a3SHansol Suh } 11646285c0a3SHansol Suh if (Hpre) { 11656285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 11666285c0a3SHansol Suh PetscCheckSameComm(tao,1,Hpre,3); 11676285c0a3SHansol Suh } 11686285c0a3SHansol Suh if (ctx) { 11696285c0a3SHansol Suh am->misfithessP = ctx; 11706285c0a3SHansol Suh } 11716285c0a3SHansol Suh if (func) { 11726285c0a3SHansol Suh am->ops->misfithess = func; 11736285c0a3SHansol Suh } 11746285c0a3SHansol Suh if (H) { 11759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H)); 11769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hx)); 11776285c0a3SHansol Suh am->Hx = H; 11786285c0a3SHansol Suh } 11796285c0a3SHansol Suh if (Hpre) { 11809566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Hpre)); 11819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hxpre)); 11826285c0a3SHansol Suh am->Hxpre = Hpre; 11836285c0a3SHansol Suh } 11846285c0a3SHansol Suh PetscFunctionReturn(0); 11856285c0a3SHansol Suh } 11866285c0a3SHansol Suh 11876285c0a3SHansol Suh /*@C 11887f5c9be9SBarry Smith TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function 11897f5c9be9SBarry Smith 11907f5c9be9SBarry Smith Collective on tao 11916285c0a3SHansol Suh 11926285c0a3SHansol Suh Input Parameters: 11936285c0a3SHansol Suh + tao - the Tao context 11946285c0a3SHansol Suh . func - function pointer for the regularizer value and gradient evaluation 11956285c0a3SHansol Suh - ctx - user context for the regularizer 11966285c0a3SHansol Suh 11976285c0a3SHansol Suh Level: advanced 1198a5a2f7acSBarry Smith 1199db781477SPatrick Sanan .seealso: `TAOADMM` 1200a5a2f7acSBarry Smith 12016285c0a3SHansol Suh @*/ 12026285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 12036285c0a3SHansol Suh { 12046285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 12056285c0a3SHansol Suh 12066285c0a3SHansol Suh PetscFunctionBegin; 12076285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 12086285c0a3SHansol Suh am->regobjgradP = ctx; 12096285c0a3SHansol Suh am->ops->regobjgrad = func; 12106285c0a3SHansol Suh PetscFunctionReturn(0); 12116285c0a3SHansol Suh } 12126285c0a3SHansol Suh 12136285c0a3SHansol Suh /*@C 12146285c0a3SHansol Suh TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back 12157f5c9be9SBarry Smith function, to be used for subsolverZ. 12167f5c9be9SBarry Smith 12177f5c9be9SBarry Smith Collective on tao 12186285c0a3SHansol Suh 12196285c0a3SHansol Suh Input Parameters: 12206285c0a3SHansol Suh + tao - the Tao context 12216285c0a3SHansol Suh . H - user-created matrix for the Hessian of the regularization term 12226285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term 12236285c0a3SHansol Suh . func - function pointer for the regularizer Hessian evaluation 12246285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 12256285c0a3SHansol Suh 12266285c0a3SHansol Suh Level: advanced 1227a5a2f7acSBarry Smith 1228db781477SPatrick Sanan .seealso: `TAOADMM` 1229a5a2f7acSBarry Smith 12306285c0a3SHansol Suh @*/ 12316285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 12326285c0a3SHansol Suh { 12336285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 12346285c0a3SHansol Suh 12356285c0a3SHansol Suh PetscFunctionBegin; 12366285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 12376285c0a3SHansol Suh if (H) { 12386285c0a3SHansol Suh PetscValidHeaderSpecific(H,MAT_CLASSID,2); 12396285c0a3SHansol Suh PetscCheckSameComm(tao,1,H,2); 12406285c0a3SHansol Suh } 12416285c0a3SHansol Suh if (Hpre) { 12426285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 12436285c0a3SHansol Suh PetscCheckSameComm(tao,1,Hpre,3); 12446285c0a3SHansol Suh } 12456285c0a3SHansol Suh if (ctx) { 12466285c0a3SHansol Suh am->reghessP = ctx; 12476285c0a3SHansol Suh } 12486285c0a3SHansol Suh if (func) { 12496285c0a3SHansol Suh am->ops->reghess = func; 12506285c0a3SHansol Suh } 12516285c0a3SHansol Suh if (H) { 12529566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H)); 12539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hz)); 12546285c0a3SHansol Suh am->Hz = H; 12556285c0a3SHansol Suh } 12566285c0a3SHansol Suh if (Hpre) { 12579566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Hpre)); 12589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hzpre)); 12596285c0a3SHansol Suh am->Hzpre = Hpre; 12606285c0a3SHansol Suh } 12616285c0a3SHansol Suh PetscFunctionReturn(0); 12626285c0a3SHansol Suh } 12636285c0a3SHansol Suh 12646285c0a3SHansol Suh /*@ 12657f5c9be9SBarry Smith TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver. 12666285c0a3SHansol Suh 12677f5c9be9SBarry Smith Collective on tao 12687f5c9be9SBarry Smith 12697f5c9be9SBarry Smith Input Parameter: 12707f5c9be9SBarry Smith . tao - the Tao context 12717f5c9be9SBarry Smith 12727f5c9be9SBarry Smith Output Parameter: 12737f5c9be9SBarry Smith . admm_tao - the parent Tao context 12746285c0a3SHansol Suh 12756285c0a3SHansol Suh Level: advanced 1276a5a2f7acSBarry Smith 1277db781477SPatrick Sanan .seealso: `TAOADMM` 1278a5a2f7acSBarry Smith 12796285c0a3SHansol Suh @*/ 12806285c0a3SHansol Suh PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao) 12816285c0a3SHansol Suh { 12826285c0a3SHansol Suh PetscFunctionBegin; 12836285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 12849566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)tao,"TaoGetADMMParentTao_ADMM", (PetscObject*) admm_tao)); 12856285c0a3SHansol Suh PetscFunctionReturn(0); 12866285c0a3SHansol Suh } 12876285c0a3SHansol Suh 12886285c0a3SHansol Suh /*@ 12897f5c9be9SBarry Smith TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state 12906285c0a3SHansol Suh 12916285c0a3SHansol Suh Not Collective 12926285c0a3SHansol Suh 12936285c0a3SHansol Suh Input Parameter: 12947f5c9be9SBarry Smith . tao - the Tao context 12957f5c9be9SBarry Smith 12967f5c9be9SBarry Smith Output Parameter: 12977f5c9be9SBarry Smith . Y - the current solution 12986285c0a3SHansol Suh 12996285c0a3SHansol Suh Level: intermediate 1300a5a2f7acSBarry Smith 1301db781477SPatrick Sanan .seealso: `TAOADMM` 1302a5a2f7acSBarry Smith 13036285c0a3SHansol Suh @*/ 13046285c0a3SHansol Suh PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y) 13056285c0a3SHansol Suh { 13066285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 13076285c0a3SHansol Suh 13086285c0a3SHansol Suh PetscFunctionBegin; 13096285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13106285c0a3SHansol Suh *Y = am->y; 13116285c0a3SHansol Suh PetscFunctionReturn(0); 13126285c0a3SHansol Suh } 13136285c0a3SHansol Suh 13146285c0a3SHansol Suh /*@ 13156285c0a3SHansol Suh TaoADMMSetRegularizerType - Set regularizer type for ADMM routine 13166285c0a3SHansol Suh 13176285c0a3SHansol Suh Not Collective 13186285c0a3SHansol Suh 1319d8d19677SJose E. Roman Input Parameters: 13206285c0a3SHansol Suh + tao - the Tao context 13217f5c9be9SBarry Smith - type - regularizer type 13227f5c9be9SBarry Smith 13237f5c9be9SBarry Smith Options Database: 1324147403d9SBarry Smith . -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer 13256285c0a3SHansol Suh 13266285c0a3SHansol Suh Level: intermediate 13276285c0a3SHansol Suh 1328db781477SPatrick Sanan .seealso: `TaoADMMGetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM` 13296285c0a3SHansol Suh @*/ 13306285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type) 13316285c0a3SHansol Suh { 13326285c0a3SHansol Suh PetscFunctionBegin; 13336285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13346285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao,type,2); 1335cac4c232SBarry Smith PetscTryMethod(tao,"TaoADMMSetRegularizerType_C",(Tao,TaoADMMRegularizerType),(tao,type)); 13366285c0a3SHansol Suh PetscFunctionReturn(0); 13376285c0a3SHansol Suh } 13386285c0a3SHansol Suh 13396285c0a3SHansol Suh /*@ 13406285c0a3SHansol Suh TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM 13416285c0a3SHansol Suh 13426285c0a3SHansol Suh Not Collective 13436285c0a3SHansol Suh 13446285c0a3SHansol Suh Input Parameter: 13456285c0a3SHansol Suh . tao - the Tao context 13466285c0a3SHansol Suh 13476285c0a3SHansol Suh Output Parameter: 13486285c0a3SHansol Suh . type - the type of regularizer 13496285c0a3SHansol Suh 13506285c0a3SHansol Suh Level: intermediate 13516285c0a3SHansol Suh 1352db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM` 13536285c0a3SHansol Suh @*/ 13546285c0a3SHansol Suh PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type) 13556285c0a3SHansol Suh { 13566285c0a3SHansol Suh PetscFunctionBegin; 13576285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1358cac4c232SBarry Smith PetscUseMethod(tao,"TaoADMMGetRegularizerType_C",(Tao,TaoADMMRegularizerType*),(tao,type)); 13596285c0a3SHansol Suh PetscFunctionReturn(0); 13606285c0a3SHansol Suh } 13616285c0a3SHansol Suh 13626285c0a3SHansol Suh /*@ 13636285c0a3SHansol Suh TaoADMMSetUpdateType - Set update routine for ADMM routine 13646285c0a3SHansol Suh 13656285c0a3SHansol Suh Not Collective 13666285c0a3SHansol Suh 1367d8d19677SJose E. Roman Input Parameters: 13686285c0a3SHansol Suh + tao - the Tao context 13696285c0a3SHansol Suh - type - spectral parameter update type 13706285c0a3SHansol Suh 13716285c0a3SHansol Suh Level: intermediate 13726285c0a3SHansol Suh 1373db781477SPatrick Sanan .seealso: `TaoADMMGetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM` 13746285c0a3SHansol Suh @*/ 13756285c0a3SHansol Suh PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type) 13766285c0a3SHansol Suh { 13776285c0a3SHansol Suh PetscFunctionBegin; 13786285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13796285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao,type,2); 1380cac4c232SBarry Smith PetscTryMethod(tao,"TaoADMMSetUpdateType_C",(Tao,TaoADMMUpdateType),(tao,type)); 13816285c0a3SHansol Suh PetscFunctionReturn(0); 13826285c0a3SHansol Suh } 13836285c0a3SHansol Suh 13846285c0a3SHansol Suh /*@ 13856285c0a3SHansol Suh TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM 13866285c0a3SHansol Suh 13876285c0a3SHansol Suh Not Collective 13886285c0a3SHansol Suh 13896285c0a3SHansol Suh Input Parameter: 1390f0fc11ceSJed Brown . tao - the Tao context 13916285c0a3SHansol Suh 13926285c0a3SHansol Suh Output Parameter: 13936285c0a3SHansol Suh . type - the type of spectral penalty update routine 13946285c0a3SHansol Suh 13956285c0a3SHansol Suh Level: intermediate 13966285c0a3SHansol Suh 1397db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM` 13986285c0a3SHansol Suh @*/ 13996285c0a3SHansol Suh PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type) 14006285c0a3SHansol Suh { 14016285c0a3SHansol Suh PetscFunctionBegin; 14026285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1403cac4c232SBarry Smith PetscUseMethod(tao,"TaoADMMGetUpdateType_C",(Tao,TaoADMMUpdateType*),(tao,type)); 14046285c0a3SHansol Suh PetscFunctionReturn(0); 14056285c0a3SHansol Suh } 1406