xref: /petsc/src/tao/constrained/impls/admm/admm.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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;
149371c9d4SSatish Balay static const char citation[] = "@misc{xu2017adaptive,\n"
156285c0a3SHansol  Suh                                "   title={Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation},\n"
166285c0a3SHansol  Suh                                "   author={Zheng Xu and Mario A. T. Figueiredo and Xiaoming Yuan and Christoph Studer and Tom Goldstein},\n"
176285c0a3SHansol  Suh                                "   year={2017},\n"
186285c0a3SHansol  Suh                                "   eprint={1704.02712},\n"
196285c0a3SHansol  Suh                                "   archivePrefix={arXiv},\n"
206285c0a3SHansol  Suh                                "   primaryClass={cs.CV}\n"
216285c0a3SHansol  Suh                                "}  \n";
226285c0a3SHansol  Suh 
23a82e8c82SStefano Zampini const char *const TaoADMMRegularizerTypes[] = {"REGULARIZER_USER", "REGULARIZER_SOFT_THRESH", "TaoADMMRegularizerType", "TAO_ADMM_", NULL};
24a82e8c82SStefano Zampini const char *const TaoADMMUpdateTypes[]      = {"UPDATE_BASIC", "UPDATE_ADAPTIVE", "UPDATE_ADAPTIVE_RELAXED", "TaoADMMUpdateType", "TAO_ADMM_", NULL};
25a82e8c82SStefano Zampini const char *const TaoALMMTypes[]            = {"CLASSIC", "PHR", "TaoALMMType", "TAO_ALMM_", NULL};
26a82e8c82SStefano Zampini 
279371c9d4SSatish Balay static PetscErrorCode TaoADMMToleranceUpdate(Tao tao) {
286285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
296285c0a3SHansol  Suh   PetscReal Axnorm, Bznorm, ATynorm, temp;
306285c0a3SHansol  Suh   Vec       tempJR, tempL;
316285c0a3SHansol  Suh   Tao       mis;
326285c0a3SHansol  Suh 
336285c0a3SHansol  Suh   PetscFunctionBegin;
346285c0a3SHansol  Suh   mis    = am->subsolverX;
356285c0a3SHansol  Suh   tempJR = am->workJacobianRight;
366285c0a3SHansol  Suh   tempL  = am->workLeft;
376285c0a3SHansol  Suh   /* ATy */
389566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianEquality(mis, am->y, mis->jacobian_equality, mis->jacobian_equality_pre));
399566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mis->jacobian_equality, am->y, tempJR));
409566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR, NORM_2, &ATynorm));
416285c0a3SHansol  Suh   /* dualres = mu * ||AT(Bz-Bzold)||_2 */
429566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR, -1., am->Bzold, am->Bz));
439566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mis->jacobian_equality, tempJR, tempL));
449566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempL, NORM_2, &am->dualres));
456285c0a3SHansol  Suh   am->dualres *= am->mu;
466285c0a3SHansol  Suh 
476285c0a3SHansol  Suh   /* ||Ax||_2, ||Bz||_2 */
489566063dSJacob Faibussowitsch   PetscCall(VecNorm(am->Ax, NORM_2, &Axnorm));
499566063dSJacob Faibussowitsch   PetscCall(VecNorm(am->Bz, NORM_2, &Bznorm));
506285c0a3SHansol  Suh 
516285c0a3SHansol  Suh   /* Set catol to be catol_admm *  max{||Ax||,||Bz||,||c||} *
526285c0a3SHansol  Suh    * Set gatol to be gatol_admm *  ||A^Ty|| *
536285c0a3SHansol  Suh    * while cnorm is ||r_k||_2, and gnorm is ||d_k||_2 */
546285c0a3SHansol  Suh   temp = am->catol_admm * PetscMax(Axnorm, (!am->const_norm) ? Bznorm : PetscMax(Bznorm, am->const_norm));
559566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintTolerances(tao, temp, PETSC_DEFAULT));
569566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao, am->gatol_admm * ATynorm, PETSC_DEFAULT, PETSC_DEFAULT));
576285c0a3SHansol  Suh   PetscFunctionReturn(0);
586285c0a3SHansol  Suh }
596285c0a3SHansol  Suh 
606285c0a3SHansol  Suh /* Penaly Update for Adaptive ADMM. */
619371c9d4SSatish Balay static PetscErrorCode AdaptiveADMMPenaltyUpdate(Tao tao) {
626285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
636285c0a3SHansol  Suh   PetscReal ydiff_norm, yhatdiff_norm, Axdiff_norm, Bzdiff_norm, Axyhat, Bzy, a_sd, a_mg, a_k, b_sd, b_mg, b_k;
646285c0a3SHansol  Suh   PetscBool hflag, gflag;
656285c0a3SHansol  Suh   Vec       tempJR, tempJR2;
666285c0a3SHansol  Suh 
676285c0a3SHansol  Suh   PetscFunctionBegin;
686285c0a3SHansol  Suh   tempJR  = am->workJacobianRight;
696285c0a3SHansol  Suh   tempJR2 = am->workJacobianRight2;
706285c0a3SHansol  Suh   hflag   = PETSC_FALSE;
716285c0a3SHansol  Suh   gflag   = PETSC_FALSE;
726285c0a3SHansol  Suh   a_k     = -1;
736285c0a3SHansol  Suh   b_k     = -1;
746285c0a3SHansol  Suh 
759566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR, -1., am->Axold, am->Ax));
769566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR2, -1., am->yhatold, am->yhat));
779566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR, NORM_2, &Axdiff_norm));
789566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR2, NORM_2, &yhatdiff_norm));
799566063dSJacob Faibussowitsch   PetscCall(VecDot(tempJR, tempJR2, &Axyhat));
806285c0a3SHansol  Suh 
819566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR, -1., am->Bz0, am->Bz));
829566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR2, -1., am->y, am->y0));
839566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR, NORM_2, &Bzdiff_norm));
849566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR2, NORM_2, &ydiff_norm));
859566063dSJacob Faibussowitsch   PetscCall(VecDot(tempJR, tempJR2, &Bzy));
866285c0a3SHansol  Suh 
876285c0a3SHansol  Suh   if (Axyhat > am->orthval * Axdiff_norm * yhatdiff_norm + am->mueps) {
886285c0a3SHansol  Suh     hflag = PETSC_TRUE;
896285c0a3SHansol  Suh     a_sd  = PetscSqr(yhatdiff_norm) / Axyhat; /* alphaSD */
906285c0a3SHansol  Suh     a_mg  = Axyhat / PetscSqr(Axdiff_norm);   /* alphaMG */
916285c0a3SHansol  Suh     a_k   = (a_mg / a_sd) > 0.5 ? a_mg : a_sd - 0.5 * a_mg;
926285c0a3SHansol  Suh   }
936285c0a3SHansol  Suh   if (Bzy > am->orthval * Bzdiff_norm * ydiff_norm + am->mueps) {
946285c0a3SHansol  Suh     gflag = PETSC_TRUE;
956285c0a3SHansol  Suh     b_sd  = PetscSqr(ydiff_norm) / Bzy;  /* betaSD */
966285c0a3SHansol  Suh     b_mg  = Bzy / PetscSqr(Bzdiff_norm); /* betaMG */
976285c0a3SHansol  Suh     b_k   = (b_mg / b_sd) > 0.5 ? b_mg : b_sd - 0.5 * b_mg;
986285c0a3SHansol  Suh   }
996285c0a3SHansol  Suh   am->muold = am->mu;
1006285c0a3SHansol  Suh   if (gflag && hflag) {
1016285c0a3SHansol  Suh     am->mu = PetscSqrtReal(a_k * b_k);
1026285c0a3SHansol  Suh   } else if (hflag) {
1036285c0a3SHansol  Suh     am->mu = a_k;
1046285c0a3SHansol  Suh   } else if (gflag) {
1056285c0a3SHansol  Suh     am->mu = b_k;
1066285c0a3SHansol  Suh   }
107ad540459SPierre Jolivet   if (am->mu > am->muold) am->mu = am->muold;
108ad540459SPierre Jolivet   if (am->mu < am->mumin) am->mu = am->mumin;
1096285c0a3SHansol  Suh   PetscFunctionReturn(0);
1106285c0a3SHansol  Suh }
1116285c0a3SHansol  Suh 
1129371c9d4SSatish Balay static PetscErrorCode TaoADMMSetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType type) {
1136285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1146285c0a3SHansol  Suh 
1156285c0a3SHansol  Suh   PetscFunctionBegin;
1166285c0a3SHansol  Suh   am->regswitch = type;
1176285c0a3SHansol  Suh   PetscFunctionReturn(0);
1186285c0a3SHansol  Suh }
1196285c0a3SHansol  Suh 
1209371c9d4SSatish Balay static PetscErrorCode TaoADMMGetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType *type) {
1216285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1226285c0a3SHansol  Suh 
1236285c0a3SHansol  Suh   PetscFunctionBegin;
1246285c0a3SHansol  Suh   *type = am->regswitch;
1256285c0a3SHansol  Suh   PetscFunctionReturn(0);
1266285c0a3SHansol  Suh }
1276285c0a3SHansol  Suh 
1289371c9d4SSatish Balay static PetscErrorCode TaoADMMSetUpdateType_ADMM(Tao tao, TaoADMMUpdateType type) {
1296285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1306285c0a3SHansol  Suh 
1316285c0a3SHansol  Suh   PetscFunctionBegin;
1326285c0a3SHansol  Suh   am->update = type;
1336285c0a3SHansol  Suh   PetscFunctionReturn(0);
1346285c0a3SHansol  Suh }
1356285c0a3SHansol  Suh 
1369371c9d4SSatish Balay static PetscErrorCode TaoADMMGetUpdateType_ADMM(Tao tao, TaoADMMUpdateType *type) {
1376285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1386285c0a3SHansol  Suh 
1396285c0a3SHansol  Suh   PetscFunctionBegin;
1406285c0a3SHansol  Suh   *type = am->update;
1416285c0a3SHansol  Suh   PetscFunctionReturn(0);
1426285c0a3SHansol  Suh }
1436285c0a3SHansol  Suh 
1446285c0a3SHansol  Suh /* This routine updates Jacobians with new x,z vectors,
1456285c0a3SHansol  Suh  * and then updates Ax and Bz vectors, then computes updated residual vector*/
1469371c9d4SSatish Balay static PetscErrorCode ADMMUpdateConstraintResidualVector(Tao tao, Vec x, Vec z, Vec Ax, Vec Bz, Vec residual) {
1476285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1486285c0a3SHansol  Suh   Tao       mis, reg;
1496285c0a3SHansol  Suh 
1506285c0a3SHansol  Suh   PetscFunctionBegin;
1516285c0a3SHansol  Suh   mis = am->subsolverX;
1526285c0a3SHansol  Suh   reg = am->subsolverZ;
1539566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianEquality(mis, x, mis->jacobian_equality, mis->jacobian_equality_pre));
1549566063dSJacob Faibussowitsch   PetscCall(MatMult(mis->jacobian_equality, x, Ax));
1559566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianEquality(reg, z, reg->jacobian_equality, reg->jacobian_equality_pre));
1569566063dSJacob Faibussowitsch   PetscCall(MatMult(reg->jacobian_equality, z, Bz));
1576285c0a3SHansol  Suh 
1589566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(residual, 1., Bz, Ax));
15948a46eb9SPierre Jolivet   if (am->constraint != NULL) PetscCall(VecAXPY(residual, -1., am->constraint));
1606285c0a3SHansol  Suh   PetscFunctionReturn(0);
1616285c0a3SHansol  Suh }
1626285c0a3SHansol  Suh 
1636285c0a3SHansol  Suh /* Updates Augmented Lagrangians to given routines *
1646285c0a3SHansol  Suh  * For subsolverX, routine needs to be ComputeObjectiveAndGraidnet
1656285c0a3SHansol  Suh  * Separate Objective and Gradient routines are not supported.  */
1669371c9d4SSatish Balay static PetscErrorCode SubObjGradUpdate(Tao tao, Vec x, PetscReal *f, Vec g, void *ptr) {
1676285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
1686285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
1696285c0a3SHansol  Suh   PetscReal temp, temp2;
1706285c0a3SHansol  Suh   Vec       tempJR;
1716285c0a3SHansol  Suh 
1726285c0a3SHansol  Suh   PetscFunctionBegin;
1736285c0a3SHansol  Suh   tempJR = am->workJacobianRight;
1749566063dSJacob Faibussowitsch   PetscCall(ADMMUpdateConstraintResidualVector(parent, x, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
1759566063dSJacob Faibussowitsch   PetscCall((*am->ops->misfitobjgrad)(am->subsolverX, x, f, g, am->misfitobjgradP));
1766285c0a3SHansol  Suh 
1776285c0a3SHansol  Suh   am->last_misfit_val = *f;
1786285c0a3SHansol  Suh   /* Objective  Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */
1799566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->y, &temp));
1809566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->residual, &temp2));
1816285c0a3SHansol  Suh   *f += temp + (am->mu / 2) * temp2;
1826285c0a3SHansol  Suh 
1836285c0a3SHansol  Suh   /* Gradient. Add + mu*AT(Ax+Bz-c) + yTA*/
1849566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian_equality, am->residual, tempJR));
1859566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, am->mu, tempJR));
1869566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian_equality, am->y, tempJR));
1879566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, 1., tempJR));
1886285c0a3SHansol  Suh   PetscFunctionReturn(0);
1896285c0a3SHansol  Suh }
1906285c0a3SHansol  Suh 
1916285c0a3SHansol  Suh /* Updates Augmented Lagrangians to given routines
1926285c0a3SHansol  Suh  * For subsolverZ, routine needs to be ComputeObjectiveAndGraidnet
1936285c0a3SHansol  Suh  * Separate Objective and Gradient routines are not supported.  */
1949371c9d4SSatish Balay static PetscErrorCode RegObjGradUpdate(Tao tao, Vec z, PetscReal *f, Vec g, void *ptr) {
1956285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
1966285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
1976285c0a3SHansol  Suh   PetscReal temp, temp2;
1986285c0a3SHansol  Suh   Vec       tempJR;
1996285c0a3SHansol  Suh 
2006285c0a3SHansol  Suh   PetscFunctionBegin;
2016285c0a3SHansol  Suh   tempJR = am->workJacobianRight;
2029566063dSJacob Faibussowitsch   PetscCall(ADMMUpdateConstraintResidualVector(parent, am->subsolverX->solution, z, am->Ax, am->Bz, am->residual));
2039566063dSJacob Faibussowitsch   PetscCall((*am->ops->regobjgrad)(am->subsolverZ, z, f, g, am->regobjgradP));
2046285c0a3SHansol  Suh   am->last_reg_val = *f;
2056285c0a3SHansol  Suh   /* Objective  Add  + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */
2069566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->y, &temp));
2079566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->residual, &temp2));
2086285c0a3SHansol  Suh   *f += temp + (am->mu / 2) * temp2;
2096285c0a3SHansol  Suh 
2106285c0a3SHansol  Suh   /* Gradient. Add + mu*BT(Ax+Bz-c) + yTB*/
2119566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality, am->residual, tempJR));
2129566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, am->mu, tempJR));
2139566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality, am->y, tempJR));
2149566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, 1., tempJR));
2156285c0a3SHansol  Suh   PetscFunctionReturn(0);
2166285c0a3SHansol  Suh }
2176285c0a3SHansol  Suh 
2186285c0a3SHansol  Suh /* Computes epsilon padded L1 norm lambda*sum(sqrt(x^2+eps^2)-eps */
2199371c9d4SSatish Balay static PetscErrorCode ADMML1EpsilonNorm(Tao tao, Vec x, PetscReal eps, PetscReal *norm) {
2206285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
2216285c0a3SHansol  Suh   PetscInt  N;
2226285c0a3SHansol  Suh 
2236285c0a3SHansol  Suh   PetscFunctionBegin;
2249566063dSJacob Faibussowitsch   PetscCall(VecGetSize(am->workLeft, &N));
2259566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(am->workLeft, x, x));
2269566063dSJacob Faibussowitsch   PetscCall(VecShift(am->workLeft, am->l1epsilon * am->l1epsilon));
2279566063dSJacob Faibussowitsch   PetscCall(VecSqrtAbs(am->workLeft));
2289566063dSJacob Faibussowitsch   PetscCall(VecSum(am->workLeft, norm));
2296285c0a3SHansol  Suh   *norm += N * am->l1epsilon;
2306285c0a3SHansol  Suh   *norm *= am->lambda;
2316285c0a3SHansol  Suh   PetscFunctionReturn(0);
2326285c0a3SHansol  Suh }
2336285c0a3SHansol  Suh 
2349371c9d4SSatish Balay static PetscErrorCode ADMMInternalHessianUpdate(Mat H, Mat Constraint, PetscBool Identity, void *ptr) {
2356285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)ptr;
2366285c0a3SHansol  Suh 
2376285c0a3SHansol  Suh   PetscFunctionBegin;
2386285c0a3SHansol  Suh   switch (am->update) {
2399371c9d4SSatish Balay   case (TAO_ADMM_UPDATE_BASIC): break;
2406285c0a3SHansol  Suh   case (TAO_ADMM_UPDATE_ADAPTIVE):
2416285c0a3SHansol  Suh   case (TAO_ADMM_UPDATE_ADAPTIVE_RELAXED):
2426285c0a3SHansol  Suh     if (H && (am->muold != am->mu)) {
2436285c0a3SHansol  Suh       if (!Identity) {
2449566063dSJacob Faibussowitsch         PetscCall(MatAXPY(H, am->mu - am->muold, Constraint, DIFFERENT_NONZERO_PATTERN));
2456285c0a3SHansol  Suh       } else {
2469566063dSJacob Faibussowitsch         PetscCall(MatShift(H, am->mu - am->muold));
2476285c0a3SHansol  Suh       }
2486285c0a3SHansol  Suh     }
2496285c0a3SHansol  Suh     break;
2506285c0a3SHansol  Suh   }
2516285c0a3SHansol  Suh   PetscFunctionReturn(0);
2526285c0a3SHansol  Suh }
2536285c0a3SHansol  Suh 
2546285c0a3SHansol  Suh /* Updates Hessian - adds second derivative of augmented Lagrangian
2556285c0a3SHansol  Suh  * H \gets H + \rho*ATA
2566285c0a3SHansol  Suh  * Here, \rho does not change in TAO_ADMM_UPDATE_BASIC - thus no-op
2576285c0a3SHansol  Suh  * For ADAPTAIVE,ADAPTIVE_RELAXED,
2586285c0a3SHansol  Suh  * H \gets H + (\rho-\rhoold)*ATA
2596285c0a3SHansol  Suh  * Here, we assume that A is linear constraint i.e., doesnt change.
2606285c0a3SHansol  Suh  * Thus, for both ADAPTIVE, and RELAXED, ATA matrix is pre-set (except for A=I (null case)) see TaoSetUp_ADMM */
2619371c9d4SSatish Balay static PetscErrorCode SubHessianUpdate(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) {
2626285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
2636285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
2646285c0a3SHansol  Suh 
2656285c0a3SHansol  Suh   PetscFunctionBegin;
2666285c0a3SHansol  Suh   if (am->Hxchange) {
2676285c0a3SHansol  Suh     /* Case where Hessian gets updated with respect to x vector input. */
2689566063dSJacob Faibussowitsch     PetscCall((*am->ops->misfithess)(am->subsolverX, x, H, Hpre, am->misfithessP));
2699566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian, am->ATA, am->xJI, am));
2706285c0a3SHansol  Suh   } else if (am->Hxbool) {
2716285c0a3SHansol  Suh     /* Hessian doesn't get updated. H(x) = c */
2726285c0a3SHansol  Suh     /* Update Lagrangian only only per TAO call */
2739566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian, am->ATA, am->xJI, am));
2746285c0a3SHansol  Suh     am->Hxbool = PETSC_FALSE;
2756285c0a3SHansol  Suh   }
2766285c0a3SHansol  Suh   PetscFunctionReturn(0);
2776285c0a3SHansol  Suh }
2786285c0a3SHansol  Suh 
2796285c0a3SHansol  Suh /* Same as SubHessianUpdate, except for B matrix instead of A matrix */
2809371c9d4SSatish Balay static PetscErrorCode RegHessianUpdate(Tao tao, Vec z, Mat H, Mat Hpre, void *ptr) {
2816285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
2826285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
2836285c0a3SHansol  Suh 
2846285c0a3SHansol  Suh   PetscFunctionBegin;
2856285c0a3SHansol  Suh 
2866285c0a3SHansol  Suh   if (am->Hzchange) {
2876285c0a3SHansol  Suh     /* Case where Hessian gets updated with respect to x vector input. */
2889566063dSJacob Faibussowitsch     PetscCall((*am->ops->reghess)(am->subsolverZ, z, H, Hpre, am->reghessP));
2899566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian, am->BTB, am->zJI, am));
2906285c0a3SHansol  Suh   } else if (am->Hzbool) {
2916285c0a3SHansol  Suh     /* Hessian doesn't get updated. H(x) = c */
2926285c0a3SHansol  Suh     /* Update Lagrangian only only per TAO call */
2939566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian, am->BTB, am->zJI, am));
2946285c0a3SHansol  Suh     am->Hzbool = PETSC_FALSE;
2956285c0a3SHansol  Suh   }
2966285c0a3SHansol  Suh   PetscFunctionReturn(0);
2976285c0a3SHansol  Suh }
2986285c0a3SHansol  Suh 
2996285c0a3SHansol  Suh /* Shell Matrix routine for A matrix.
3006285c0a3SHansol  Suh  * This gets used when user puts NULL for
3016285c0a3SHansol  Suh  * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...)
3026285c0a3SHansol  Suh  * Essentially sets A=I*/
3039371c9d4SSatish Balay static PetscErrorCode JacobianIdentity(Mat mat, Vec in, Vec out) {
3046285c0a3SHansol  Suh   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(VecCopy(in, out));
3066285c0a3SHansol  Suh   PetscFunctionReturn(0);
3076285c0a3SHansol  Suh }
3086285c0a3SHansol  Suh 
3096285c0a3SHansol  Suh /* Shell Matrix routine for B matrix.
3106285c0a3SHansol  Suh  * This gets used when user puts NULL for
3116285c0a3SHansol  Suh  * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...)
3126285c0a3SHansol  Suh  * Sets B=-I */
3139371c9d4SSatish Balay static PetscErrorCode JacobianIdentityB(Mat mat, Vec in, Vec out) {
3146285c0a3SHansol  Suh   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(VecCopy(in, out));
3169566063dSJacob Faibussowitsch   PetscCall(VecScale(out, -1.));
3176285c0a3SHansol  Suh   PetscFunctionReturn(0);
3186285c0a3SHansol  Suh }
3196285c0a3SHansol  Suh 
3206285c0a3SHansol  Suh /* Solve f(x) + g(z) s.t. Ax + Bz = c */
3219371c9d4SSatish Balay static PetscErrorCode TaoSolve_ADMM(Tao tao) {
3226285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
3236285c0a3SHansol  Suh   PetscInt  N;
3246285c0a3SHansol  Suh   PetscReal reg_func;
3256285c0a3SHansol  Suh   PetscBool is_reg_shell;
3266285c0a3SHansol  Suh   Vec       tempL;
3276285c0a3SHansol  Suh 
3286285c0a3SHansol  Suh   PetscFunctionBegin;
3296285c0a3SHansol  Suh   if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) {
3303c859ba3SBarry Smith     PetscCheck(am->subsolverX->ops->computejacobianequality, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoADMMSetMisfitConstraintJacobian() first");
3313c859ba3SBarry Smith     PetscCheck(am->subsolverZ->ops->computejacobianequality, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoADMMSetRegularizerConstraintJacobian() first");
33248a46eb9SPierre Jolivet     if (am->constraint != NULL) PetscCall(VecNorm(am->constraint, NORM_2, &am->const_norm));
3336285c0a3SHansol  Suh   }
3346285c0a3SHansol  Suh   tempL = am->workLeft;
3359566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tempL, &N));
3366285c0a3SHansol  Suh 
33748a46eb9SPierre Jolivet   if (am->Hx && am->ops->misfithess) PetscCall(TaoSetHessian(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao));
3386285c0a3SHansol  Suh 
3396285c0a3SHansol  Suh   if (!am->zJI) {
3406285c0a3SHansol  Suh     /* Currently, B is assumed to be a linear system, i.e., not getting updated*/
3419566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(am->JB, am->JB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(am->BTB)));
3426285c0a3SHansol  Suh   }
3436285c0a3SHansol  Suh   if (!am->xJI) {
3446285c0a3SHansol  Suh     /* Currently, A is assumed to be a linear system, i.e., not getting updated*/
3459566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(am->ATA)));
3466285c0a3SHansol  Suh   }
3476285c0a3SHansol  Suh 
3486285c0a3SHansol  Suh   is_reg_shell = PETSC_FALSE;
3496285c0a3SHansol  Suh 
3509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell));
3516285c0a3SHansol  Suh 
3526285c0a3SHansol  Suh   if (!is_reg_shell) {
3536285c0a3SHansol  Suh     switch (am->regswitch) {
3549371c9d4SSatish Balay     case (TAO_ADMM_REGULARIZER_USER): break;
3556285c0a3SHansol  Suh     case (TAO_ADMM_REGULARIZER_SOFT_THRESH):
3566285c0a3SHansol  Suh       /* Soft Threshold. */
3576285c0a3SHansol  Suh       break;
3586285c0a3SHansol  Suh     }
3591baa6e33SBarry Smith     if (am->ops->regobjgrad) PetscCall(TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao));
36048a46eb9SPierre Jolivet     if (am->Hz && am->ops->reghess) PetscCall(TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao));
3616285c0a3SHansol  Suh   }
3626285c0a3SHansol  Suh 
3636285c0a3SHansol  Suh   switch (am->update) {
3646285c0a3SHansol  Suh   case TAO_ADMM_UPDATE_BASIC:
3656285c0a3SHansol  Suh     if (am->subsolverX->hessian) {
3666285c0a3SHansol  Suh       /* In basic case, Hessian does not get updated w.r.t. to spectral penalty
3676285c0a3SHansol  Suh        * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/
3686285c0a3SHansol  Suh       if (!am->xJI) {
3699566063dSJacob Faibussowitsch         PetscCall(MatAXPY(am->subsolverX->hessian, am->mu, am->ATA, DIFFERENT_NONZERO_PATTERN));
3706285c0a3SHansol  Suh       } else {
3719566063dSJacob Faibussowitsch         PetscCall(MatShift(am->subsolverX->hessian, am->mu));
3726285c0a3SHansol  Suh       }
3736285c0a3SHansol  Suh     }
3746285c0a3SHansol  Suh     if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) {
3756285c0a3SHansol  Suh       if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) {
3769566063dSJacob Faibussowitsch         PetscCall(MatAXPY(am->subsolverZ->hessian, am->mu, am->BTB, DIFFERENT_NONZERO_PATTERN));
3776285c0a3SHansol  Suh       } else {
3789566063dSJacob Faibussowitsch         PetscCall(MatShift(am->subsolverZ->hessian, am->mu));
3796285c0a3SHansol  Suh       }
3806285c0a3SHansol  Suh     }
3816285c0a3SHansol  Suh     break;
3826285c0a3SHansol  Suh   case TAO_ADMM_UPDATE_ADAPTIVE:
3839371c9d4SSatish Balay   case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: break;
3846285c0a3SHansol  Suh   }
3856285c0a3SHansol  Suh 
3869566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
3876285c0a3SHansol  Suh   tao->reason = TAO_CONTINUE_ITERATING;
3886285c0a3SHansol  Suh 
3896285c0a3SHansol  Suh   while (tao->reason == TAO_CONTINUE_ITERATING) {
390dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
3919566063dSJacob Faibussowitsch     PetscCall(VecCopy(am->Bz, am->Bzold));
3926285c0a3SHansol  Suh 
3936285c0a3SHansol  Suh     /* x update */
3949566063dSJacob Faibussowitsch     PetscCall(TaoSolve(am->subsolverX));
3959566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre));
3969566063dSJacob Faibussowitsch     PetscCall(MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution, am->Ax));
3976285c0a3SHansol  Suh 
3986285c0a3SHansol  Suh     am->Hxbool = PETSC_TRUE;
3996285c0a3SHansol  Suh 
4006285c0a3SHansol  Suh     /* z update */
4016285c0a3SHansol  Suh     switch (am->regswitch) {
4029371c9d4SSatish Balay     case TAO_ADMM_REGULARIZER_USER: PetscCall(TaoSolve(am->subsolverZ)); break;
4036285c0a3SHansol  Suh     case TAO_ADMM_REGULARIZER_SOFT_THRESH:
4046285c0a3SHansol  Suh       /* L1 assumes A,B jacobians are identity nxn matrix */
4059566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(am->workJacobianRight, 1 / am->mu, am->y, am->Ax));
4069566063dSJacob Faibussowitsch       PetscCall(TaoSoftThreshold(am->workJacobianRight, -am->lambda / am->mu, am->lambda / am->mu, am->subsolverZ->solution));
4076285c0a3SHansol  Suh       break;
4086285c0a3SHansol  Suh     }
4096285c0a3SHansol  Suh     am->Hzbool = PETSC_TRUE;
4106285c0a3SHansol  Suh     /* Returns Ax + Bz - c with updated Ax,Bz vectors */
4119566063dSJacob Faibussowitsch     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
4126285c0a3SHansol  Suh     /* Dual variable, y += y + mu*(Ax+Bz-c) */
4139566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(am->y, am->mu, am->residual, am->yold));
4146285c0a3SHansol  Suh 
4156285c0a3SHansol  Suh     /* stopping tolerance update */
4169566063dSJacob Faibussowitsch     PetscCall(TaoADMMToleranceUpdate(tao));
4176285c0a3SHansol  Suh 
4186285c0a3SHansol  Suh     /* Updating Spectral Penalty */
4196285c0a3SHansol  Suh     switch (am->update) {
4209371c9d4SSatish Balay     case TAO_ADMM_UPDATE_BASIC: am->muold = am->mu; break;
4216285c0a3SHansol  Suh     case TAO_ADMM_UPDATE_ADAPTIVE:
4226285c0a3SHansol  Suh     case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
4236285c0a3SHansol  Suh       if (tao->niter == 0) {
4249566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->y, am->y0));
4259566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
4261baa6e33SBarry Smith         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
4279566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold));
4289566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Ax, am->Axold));
4299566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Bz, am->Bz0));
4306285c0a3SHansol  Suh         am->muold = am->mu;
4316285c0a3SHansol  Suh       } else if (tao->niter % am->T == 1) {
4326285c0a3SHansol  Suh         /* we have compute Bzold in a previous iteration, and we computed Ax above */
4339566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
4341baa6e33SBarry Smith         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
4359566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->yhat, -am->mu, am->residual, am->yold));
4369566063dSJacob Faibussowitsch         PetscCall(AdaptiveADMMPenaltyUpdate(tao));
4379566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Ax, am->Axold));
4389566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Bz, am->Bz0));
4399566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->yhat, am->yhatold));
4409566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->y, am->y0));
4416285c0a3SHansol  Suh       } else {
4426285c0a3SHansol  Suh         am->muold = am->mu;
4436285c0a3SHansol  Suh       }
4446285c0a3SHansol  Suh       break;
4459371c9d4SSatish Balay     default: break;
4466285c0a3SHansol  Suh     }
4476285c0a3SHansol  Suh     tao->niter++;
4486285c0a3SHansol  Suh 
4496285c0a3SHansol  Suh     /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/
4506285c0a3SHansol  Suh     switch (am->regswitch) {
4516285c0a3SHansol  Suh     case TAO_ADMM_REGULARIZER_USER:
4526285c0a3SHansol  Suh       if (is_reg_shell) {
4539566063dSJacob Faibussowitsch         PetscCall(ADMML1EpsilonNorm(tao, am->subsolverZ->solution, am->l1epsilon, &reg_func));
4546285c0a3SHansol  Suh       } else {
4556285c0a3SHansol  Suh         (*am->ops->regobjgrad)(am->subsolverZ, am->subsolverX->solution, &reg_func, tempL, am->regobjgradP);
4566285c0a3SHansol  Suh       }
4576285c0a3SHansol  Suh       break;
4589371c9d4SSatish Balay     case TAO_ADMM_REGULARIZER_SOFT_THRESH: PetscCall(ADMML1EpsilonNorm(tao, am->subsolverZ->solution, am->l1epsilon, &reg_func)); break;
4596285c0a3SHansol  Suh     }
4609566063dSJacob Faibussowitsch     PetscCall(VecCopy(am->y, am->yold));
4619566063dSJacob Faibussowitsch     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
4629566063dSJacob Faibussowitsch     PetscCall(VecNorm(am->residual, NORM_2, &am->resnorm));
4639566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, am->last_misfit_val + reg_func, am->dualres, am->resnorm, tao->ksp_its));
4646285c0a3SHansol  Suh 
4659566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, am->last_misfit_val + reg_func, am->dualres, am->resnorm, 1.0));
466dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
4676285c0a3SHansol  Suh   }
4686285c0a3SHansol  Suh   /* Update vectors */
4699566063dSJacob Faibussowitsch   PetscCall(VecCopy(am->subsolverX->solution, tao->solution));
4709566063dSJacob Faibussowitsch   PetscCall(VecCopy(am->subsolverX->gradient, tao->gradient));
4719566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX, "TaoGetADMMParentTao_ADMM", NULL));
4729566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ, "TaoGetADMMParentTao_ADMM", NULL));
4739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", NULL));
4749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", NULL));
4759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", NULL));
4769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", NULL));
4776285c0a3SHansol  Suh   PetscFunctionReturn(0);
4786285c0a3SHansol  Suh }
4796285c0a3SHansol  Suh 
4809371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_ADMM(Tao tao, PetscOptionItems *PetscOptionsObject) {
4816285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
4826285c0a3SHansol  Suh 
4836285c0a3SHansol  Suh   PetscFunctionBegin;
484d0609cedSBarry 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. ");
4859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient", "regularizer constant", "", am->lambda, &am->lambda, NULL));
4869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty", "Constant for Augmented Lagrangian term.", "", am->mu, &am->mu, NULL));
4879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter", "x relaxation parameter for Z update.", "", am->gamma, &am->gamma, NULL));
4889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor", "ADMM dynamic tolerance update factor.", "", am->tol, &am->tol, NULL));
4899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor", "ADMM spectral penalty update curvature safeguard value.", "", am->orthval, &am->orthval, NULL));
4909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty", "Set ADMM minimum spectral penalty.", "", am->mumin, &am->mumin, NULL));
491d0609cedSBarry Smith   PetscCall(PetscOptionsEnum("-tao_admm_dual_update", "Lagrangian dual update policy", "TaoADMMUpdateType", TaoADMMUpdateTypes, (PetscEnum)am->update, (PetscEnum *)&am->update, NULL));
492d0609cedSBarry Smith   PetscCall(PetscOptionsEnum("-tao_admm_regularizer_type", "ADMM regularizer update rule", "TaoADMMRegularizerType", TaoADMMRegularizerTypes, (PetscEnum)am->regswitch, (PetscEnum *)&am->regswitch, NULL));
493d0609cedSBarry Smith   PetscOptionsHeadEnd();
4949566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(am->subsolverX));
49548a46eb9SPierre Jolivet   if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) PetscCall(TaoSetFromOptions(am->subsolverZ));
4966285c0a3SHansol  Suh   PetscFunctionReturn(0);
4976285c0a3SHansol  Suh }
4986285c0a3SHansol  Suh 
4999371c9d4SSatish Balay static PetscErrorCode TaoView_ADMM(Tao tao, PetscViewer viewer) {
5006285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
5016285c0a3SHansol  Suh 
5026285c0a3SHansol  Suh   PetscFunctionBegin;
5039566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
5049566063dSJacob Faibussowitsch   PetscCall(TaoView(am->subsolverX, viewer));
5059566063dSJacob Faibussowitsch   PetscCall(TaoView(am->subsolverZ, viewer));
5069566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
5076285c0a3SHansol  Suh   PetscFunctionReturn(0);
5086285c0a3SHansol  Suh }
5096285c0a3SHansol  Suh 
5109371c9d4SSatish Balay static PetscErrorCode TaoSetUp_ADMM(Tao tao) {
5116285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
5126285c0a3SHansol  Suh   PetscInt  n, N, M;
5136285c0a3SHansol  Suh 
5146285c0a3SHansol  Suh   PetscFunctionBegin;
5159566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution, &n));
5169566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &N));
5176285c0a3SHansol  Suh   /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */
5186285c0a3SHansol  Suh   if (!am->JB) {
5196285c0a3SHansol  Suh     am->zJI = PETSC_TRUE;
5209566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao), n, n, PETSC_DETERMINE, PETSC_DETERMINE, NULL, &am->JB));
5219566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JB, MATOP_MULT, (void (*)(void))JacobianIdentityB));
5229566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JB, MATOP_MULT_TRANSPOSE, (void (*)(void))JacobianIdentityB));
5236285c0a3SHansol  Suh     am->JBpre = am->JB;
5246285c0a3SHansol  Suh   }
5256285c0a3SHansol  Suh   if (!am->JA) {
5266285c0a3SHansol  Suh     am->xJI = PETSC_TRUE;
5279566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao), n, n, PETSC_DETERMINE, PETSC_DETERMINE, NULL, &am->JA));
5289566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JA, MATOP_MULT, (void (*)(void))JacobianIdentity));
5299566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JA, MATOP_MULT_TRANSPOSE, (void (*)(void))JacobianIdentity));
5306285c0a3SHansol  Suh     am->JApre = am->JA;
5316285c0a3SHansol  Suh   }
5329566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(am->JA, NULL, &am->Ax));
53348a46eb9SPierre Jolivet   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
5349566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(am->subsolverX, tao->solution));
5356285c0a3SHansol  Suh   if (!am->z) {
5369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &am->z));
5379566063dSJacob Faibussowitsch     PetscCall(VecSet(am->z, 0.0));
5386285c0a3SHansol  Suh   }
5399566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(am->subsolverZ, am->z));
54048a46eb9SPierre Jolivet   if (!am->workLeft) PetscCall(VecDuplicate(tao->solution, &am->workLeft));
54148a46eb9SPierre Jolivet   if (!am->Axold) PetscCall(VecDuplicate(am->Ax, &am->Axold));
54248a46eb9SPierre Jolivet   if (!am->workJacobianRight) PetscCall(VecDuplicate(am->Ax, &am->workJacobianRight));
54348a46eb9SPierre Jolivet   if (!am->workJacobianRight2) PetscCall(VecDuplicate(am->Ax, &am->workJacobianRight2));
54448a46eb9SPierre Jolivet   if (!am->Bz) PetscCall(VecDuplicate(am->Ax, &am->Bz));
54548a46eb9SPierre Jolivet   if (!am->Bzold) PetscCall(VecDuplicate(am->Ax, &am->Bzold));
54648a46eb9SPierre Jolivet   if (!am->Bz0) PetscCall(VecDuplicate(am->Ax, &am->Bz0));
5476285c0a3SHansol  Suh   if (!am->y) {
5489566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->y));
5499566063dSJacob Faibussowitsch     PetscCall(VecSet(am->y, 0.0));
5506285c0a3SHansol  Suh   }
5516285c0a3SHansol  Suh   if (!am->yold) {
5529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->yold));
5539566063dSJacob Faibussowitsch     PetscCall(VecSet(am->yold, 0.0));
5546285c0a3SHansol  Suh   }
5556285c0a3SHansol  Suh   if (!am->y0) {
5569566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->y0));
5579566063dSJacob Faibussowitsch     PetscCall(VecSet(am->y0, 0.0));
5586285c0a3SHansol  Suh   }
5596285c0a3SHansol  Suh   if (!am->yhat) {
5609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->yhat));
5619566063dSJacob Faibussowitsch     PetscCall(VecSet(am->yhat, 0.0));
5626285c0a3SHansol  Suh   }
5636285c0a3SHansol  Suh   if (!am->yhatold) {
5649566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->yhatold));
5659566063dSJacob Faibussowitsch     PetscCall(VecSet(am->yhatold, 0.0));
5666285c0a3SHansol  Suh   }
5676285c0a3SHansol  Suh   if (!am->residual) {
5689566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->residual));
5699566063dSJacob Faibussowitsch     PetscCall(VecSet(am->residual, 0.0));
5706285c0a3SHansol  Suh   }
5716285c0a3SHansol  Suh   if (!am->constraint) {
5726285c0a3SHansol  Suh     am->constraint = NULL;
5736285c0a3SHansol  Suh   } else {
5749566063dSJacob Faibussowitsch     PetscCall(VecGetSize(am->constraint, &M));
5753c859ba3SBarry Smith     PetscCheck(M == N, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Solution vector and constraint vector must be of same size!");
5766285c0a3SHansol  Suh   }
5776285c0a3SHansol  Suh 
5786285c0a3SHansol  Suh   /* Save changed tao tolerance for adaptive tolerance */
579ad540459SPierre Jolivet   if (tao->gatol_changed) am->gatol_admm = tao->gatol;
580ad540459SPierre Jolivet   if (tao->catol_changed) am->catol_admm = tao->catol;
5816285c0a3SHansol  Suh 
5826285c0a3SHansol  Suh   /*Update spectral and dual elements to X subsolver */
5839566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao));
5849566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX, am->JA, am->JApre, am->ops->misfitjac, am->misfitjacobianP));
5859566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ, am->JB, am->JBpre, am->ops->regjac, am->regjacobianP));
5866285c0a3SHansol  Suh   PetscFunctionReturn(0);
5876285c0a3SHansol  Suh }
5886285c0a3SHansol  Suh 
5899371c9d4SSatish Balay static PetscErrorCode TaoDestroy_ADMM(Tao tao) {
5906285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
5916285c0a3SHansol  Suh 
5926285c0a3SHansol  Suh   PetscFunctionBegin;
5939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->z));
5949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Ax));
5959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Axold));
5969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Bz));
5979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Bzold));
5989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Bz0));
5999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->residual));
6009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->y));
6019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->yold));
6029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->y0));
6039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->yhat));
6049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->yhatold));
6059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->workLeft));
6069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->workJacobianRight));
6079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->workJacobianRight2));
6086285c0a3SHansol  Suh 
6099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->JA));
6109566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->JB));
61148a46eb9SPierre Jolivet   if (!am->xJI) PetscCall(MatDestroy(&am->JApre));
61248a46eb9SPierre Jolivet   if (!am->zJI) PetscCall(MatDestroy(&am->JBpre));
6136285c0a3SHansol  Suh   if (am->Hx) {
6149566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hx));
6159566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hxpre));
6166285c0a3SHansol  Suh   }
6176285c0a3SHansol  Suh   if (am->Hz) {
6189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hz));
6199566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hzpre));
6206285c0a3SHansol  Suh   }
6219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->ATA));
6229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->BTB));
6239566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&am->subsolverX));
6249566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&am->subsolverZ));
6256285c0a3SHansol  Suh   am->parent = NULL;
6262e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", NULL));
6272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", NULL));
6282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", NULL));
6292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", NULL));
6309566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
6316285c0a3SHansol  Suh   PetscFunctionReturn(0);
6326285c0a3SHansol  Suh }
6336285c0a3SHansol  Suh 
6346285c0a3SHansol  Suh /*MC
6356285c0a3SHansol  Suh 
6366285c0a3SHansol  Suh   TAOADMM - Alternating direction method of multipliers method fo solving linear problems with
6376285c0a3SHansol  Suh             constraints. in a min_x f(x) + g(z)  s.t. Ax+Bz=c.
638a5b23f4aSJose E. Roman             This algorithm employs two sub Tao solvers, of which type can be specified
6396285c0a3SHansol  Suh             by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers.
6406285c0a3SHansol  Suh             Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via
6416285c0a3SHansol  Suh             TaoADMMSet{Misfit,Regularizer}HessianChangeStatus.
6426285c0a3SHansol  Suh             Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type.
6436285c0a3SHansol  Suh             There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update,
644a5b23f4aSJose E. Roman             currently there are basic option and adaptive option.
6456285c0a3SHansol  Suh             Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian.
6466285c0a3SHansol  Suh             c can be set with TaoADMMSetConstraintVectorRHS.
6476285c0a3SHansol  Suh             The user can also provide regularizer weight for second subsolver.
6486285c0a3SHansol  Suh 
6496285c0a3SHansol  Suh   References:
650606c0280SSatish Balay . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom
6516285c0a3SHansol  Suh           "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation"
6526285c0a3SHansol  Suh           The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017.
6536285c0a3SHansol  Suh 
6546285c0a3SHansol  Suh   Options Database Keys:
6556285c0a3SHansol  Suh + -tao_admm_regularizer_coefficient        - regularizer constant (default 1.e-6)
6566285c0a3SHansol  Suh . -tao_admm_spectral_penalty               - Constant for Augmented Lagrangian term (default 1.)
6576285c0a3SHansol  Suh . -tao_admm_relaxation_parameter           - relaxation parameter for Z update (default 1.)
6586285c0a3SHansol  Suh . -tao_admm_tolerance_update_factor        - ADMM dynamic tolerance update factor (default 1.e-12)
6596285c0a3SHansol  Suh . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2)
6606285c0a3SHansol  Suh . -tao_admm_minimum_spectral_penalty       - Set ADMM minimum spectral penalty (default 0)
6616285c0a3SHansol  Suh . -tao_admm_dual_update                    - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic")
6626285c0a3SHansol  Suh - -tao_admm_regularizer_type               - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold")
6636285c0a3SHansol  Suh 
6646285c0a3SHansol  Suh   Level: beginner
6656285c0a3SHansol  Suh 
666db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitHessianChangeStatus()`, `TaoADMMSetRegHessianChangeStatus()`, `TaoADMMGetSpectralPenalty()`,
667db781477SPatrick Sanan           `TaoADMMGetMisfitSubsolver()`, `TaoADMMGetRegularizationSubsolver()`, `TaoADMMSetConstraintVectorRHS()`,
668db781477SPatrick Sanan           `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetRegularizerCoefficient()`,
669db781477SPatrick Sanan           `TaoADMMSetRegularizerConstraintJacobian()`, `TaoADMMSetMisfitConstraintJacobian()`,
670db781477SPatrick Sanan           `TaoADMMSetMisfitObjectiveAndGradientRoutine()`, `TaoADMMSetMisfitHessianRoutine()`,
671db781477SPatrick Sanan           `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, `TaoADMMSetRegularizerHessianRoutine()`,
672db781477SPatrick Sanan           `TaoGetADMMParentTao()`, `TaoADMMGetDualVector()`, `TaoADMMSetRegularizerType()`,
673db781477SPatrick Sanan           `TaoADMMGetRegularizerType()`, `TaoADMMSetUpdateType()`, `TaoADMMGetUpdateType()`
6746285c0a3SHansol  Suh M*/
6756285c0a3SHansol  Suh 
6769371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao) {
6776285c0a3SHansol  Suh   TAO_ADMM *am;
6786285c0a3SHansol  Suh 
6796285c0a3SHansol  Suh   PetscFunctionBegin;
680*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&am));
6816285c0a3SHansol  Suh 
6826285c0a3SHansol  Suh   tao->ops->destroy        = TaoDestroy_ADMM;
6836285c0a3SHansol  Suh   tao->ops->setup          = TaoSetUp_ADMM;
6846285c0a3SHansol  Suh   tao->ops->setfromoptions = TaoSetFromOptions_ADMM;
6856285c0a3SHansol  Suh   tao->ops->view           = TaoView_ADMM;
6866285c0a3SHansol  Suh   tao->ops->solve          = TaoSolve_ADMM;
6876285c0a3SHansol  Suh 
6886285c0a3SHansol  Suh   tao->data           = (void *)am;
6896285c0a3SHansol  Suh   am->l1epsilon       = 1e-6;
6906285c0a3SHansol  Suh   am->lambda          = 1e-4;
6916285c0a3SHansol  Suh   am->mu              = 1.;
6926285c0a3SHansol  Suh   am->muold           = 0.;
6936285c0a3SHansol  Suh   am->mueps           = PETSC_MACHINE_EPSILON;
6946285c0a3SHansol  Suh   am->mumin           = 0.;
6956285c0a3SHansol  Suh   am->orthval         = 0.2;
6966285c0a3SHansol  Suh   am->T               = 2;
6976285c0a3SHansol  Suh   am->parent          = tao;
6986285c0a3SHansol  Suh   am->update          = TAO_ADMM_UPDATE_BASIC;
6996285c0a3SHansol  Suh   am->regswitch       = TAO_ADMM_REGULARIZER_SOFT_THRESH;
7006285c0a3SHansol  Suh   am->tol             = PETSC_SMALL;
7016285c0a3SHansol  Suh   am->const_norm      = 0;
7026285c0a3SHansol  Suh   am->resnorm         = 0;
7036285c0a3SHansol  Suh   am->dualres         = 0;
70483c8fe1dSLisandro Dalcin   am->ops->regobjgrad = NULL;
70583c8fe1dSLisandro Dalcin   am->ops->reghess    = NULL;
7066285c0a3SHansol  Suh   am->gamma           = 1;
70783c8fe1dSLisandro Dalcin   am->regobjgradP     = NULL;
70883c8fe1dSLisandro Dalcin   am->reghessP        = NULL;
7096285c0a3SHansol  Suh   am->gatol_admm      = 1e-8;
7106285c0a3SHansol  Suh   am->catol_admm      = 0;
7116285c0a3SHansol  Suh   am->Hxchange        = PETSC_TRUE;
7126285c0a3SHansol  Suh   am->Hzchange        = PETSC_TRUE;
7136285c0a3SHansol  Suh   am->Hzbool          = PETSC_TRUE;
7146285c0a3SHansol  Suh   am->Hxbool          = PETSC_TRUE;
7156285c0a3SHansol  Suh 
7169566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &am->subsolverX));
7179566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(am->subsolverX, "misfit_"));
7189566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX, (PetscObject)tao, 1));
7199566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &am->subsolverZ));
7209566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(am->subsolverZ, "reg_"));
7219566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ, (PetscObject)tao, 1));
7226285c0a3SHansol  Suh 
7239566063dSJacob Faibussowitsch   PetscCall(TaoSetType(am->subsolverX, TAONLS));
7249566063dSJacob Faibussowitsch   PetscCall(TaoSetType(am->subsolverZ, TAONLS));
7259566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX, "TaoGetADMMParentTao_ADMM", (PetscObject)tao));
7269566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ, "TaoGetADMMParentTao_ADMM", (PetscObject)tao));
7279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", TaoADMMSetRegularizerType_ADMM));
7289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", TaoADMMGetRegularizerType_ADMM));
7299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", TaoADMMSetUpdateType_ADMM));
7309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", TaoADMMGetUpdateType_ADMM));
7316285c0a3SHansol  Suh   PetscFunctionReturn(0);
7326285c0a3SHansol  Suh }
7336285c0a3SHansol  Suh 
7346285c0a3SHansol  Suh /*@
7356285c0a3SHansol  Suh   TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines  whether Hessian matrix of misfit subsolver changes with respect to input vector.
7366285c0a3SHansol  Suh 
7376285c0a3SHansol  Suh   Collective on Tao
7386285c0a3SHansol  Suh 
7396285c0a3SHansol  Suh   Input Parameters:
7407f5c9be9SBarry Smith +  tao - the Tao solver context.
7417f5c9be9SBarry Smith -  b - the Hessian matrix change status boolean, PETSC_FALSE  when the Hessian matrix does not change, TRUE otherwise.
7426285c0a3SHansol  Suh 
7436285c0a3SHansol  Suh   Level: advanced
744a5a2f7acSBarry Smith 
745db781477SPatrick Sanan .seealso: `TAOADMM`
746a5a2f7acSBarry Smith 
7476285c0a3SHansol  Suh @*/
7489371c9d4SSatish Balay PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b) {
7496285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
7506285c0a3SHansol  Suh 
7516285c0a3SHansol  Suh   PetscFunctionBegin;
7526285c0a3SHansol  Suh   am->Hxchange = b;
7536285c0a3SHansol  Suh   PetscFunctionReturn(0);
7546285c0a3SHansol  Suh }
7556285c0a3SHansol  Suh 
7566285c0a3SHansol  Suh /*@
7576285c0a3SHansol  Suh   TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector.
7587f5c9be9SBarry Smith 
7596285c0a3SHansol  Suh   Collective on Tao
7606285c0a3SHansol  Suh 
7616285c0a3SHansol  Suh   Input Parameters:
7626285c0a3SHansol  Suh +  tao - the Tao solver context
7637f5c9be9SBarry Smith -  b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise.
7646285c0a3SHansol  Suh 
7656285c0a3SHansol  Suh   Level: advanced
766a5a2f7acSBarry Smith 
767db781477SPatrick Sanan .seealso: `TAOADMM`
768a5a2f7acSBarry Smith 
7696285c0a3SHansol  Suh @*/
7709371c9d4SSatish Balay PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b) {
7716285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
7726285c0a3SHansol  Suh 
7736285c0a3SHansol  Suh   PetscFunctionBegin;
7746285c0a3SHansol  Suh   am->Hzchange = b;
7756285c0a3SHansol  Suh   PetscFunctionReturn(0);
7766285c0a3SHansol  Suh }
7776285c0a3SHansol  Suh 
7786285c0a3SHansol  Suh /*@
7796285c0a3SHansol  Suh   TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value
7806285c0a3SHansol  Suh 
7816285c0a3SHansol  Suh   Collective on Tao
7826285c0a3SHansol  Suh 
7836285c0a3SHansol  Suh   Input Parameters:
7846285c0a3SHansol  Suh +  tao - the Tao solver context
7857f5c9be9SBarry Smith -  mu - spectral penalty
7866285c0a3SHansol  Suh 
7876285c0a3SHansol  Suh   Level: advanced
7886285c0a3SHansol  Suh 
789db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TAOADMM`
7906285c0a3SHansol  Suh @*/
7919371c9d4SSatish Balay PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu) {
7926285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
7936285c0a3SHansol  Suh 
7946285c0a3SHansol  Suh   PetscFunctionBegin;
7956285c0a3SHansol  Suh   am->mu = mu;
7966285c0a3SHansol  Suh   PetscFunctionReturn(0);
7976285c0a3SHansol  Suh }
7986285c0a3SHansol  Suh 
7996285c0a3SHansol  Suh /*@
8006285c0a3SHansol  Suh   TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value
8016285c0a3SHansol  Suh 
8026285c0a3SHansol  Suh   Collective on Tao
8036285c0a3SHansol  Suh 
8047f5c9be9SBarry Smith   Input Parameter:
8057f5c9be9SBarry Smith .  tao - the Tao solver context
8067f5c9be9SBarry Smith 
8077f5c9be9SBarry Smith   Output Parameter:
8087f5c9be9SBarry Smith .  mu - spectral penalty
8096285c0a3SHansol  Suh 
8106285c0a3SHansol  Suh   Level: advanced
8116285c0a3SHansol  Suh 
812db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetSpectralPenalty()`, `TAOADMM`
8136285c0a3SHansol  Suh @*/
8149371c9d4SSatish Balay PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu) {
8156285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8166285c0a3SHansol  Suh 
8176285c0a3SHansol  Suh   PetscFunctionBegin;
8186285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
8196285c0a3SHansol  Suh   PetscValidRealPointer(mu, 2);
8206285c0a3SHansol  Suh   *mu = am->mu;
8216285c0a3SHansol  Suh   PetscFunctionReturn(0);
8226285c0a3SHansol  Suh }
8236285c0a3SHansol  Suh 
8246285c0a3SHansol  Suh /*@
8256285c0a3SHansol  Suh   TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM
8266285c0a3SHansol  Suh 
8276285c0a3SHansol  Suh   Collective on Tao
8286285c0a3SHansol  Suh 
8297f5c9be9SBarry Smith   Input Parameter:
8307f5c9be9SBarry Smith .  tao - the Tao solver context
8317f5c9be9SBarry Smith 
8327f5c9be9SBarry Smith    Output Parameter:
8337f5c9be9SBarry Smith .  misfit - the Tao subsolver context
8346285c0a3SHansol  Suh 
8356285c0a3SHansol  Suh   Level: advanced
836a5a2f7acSBarry Smith 
837db781477SPatrick Sanan .seealso: `TAOADMM`
838a5a2f7acSBarry Smith 
8396285c0a3SHansol  Suh @*/
8409371c9d4SSatish Balay PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit) {
8416285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8426285c0a3SHansol  Suh 
8436285c0a3SHansol  Suh   PetscFunctionBegin;
8446285c0a3SHansol  Suh   *misfit = am->subsolverX;
8456285c0a3SHansol  Suh   PetscFunctionReturn(0);
8466285c0a3SHansol  Suh }
8476285c0a3SHansol  Suh 
8486285c0a3SHansol  Suh /*@
8496285c0a3SHansol  Suh   TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM
8506285c0a3SHansol  Suh 
8516285c0a3SHansol  Suh   Collective on Tao
8526285c0a3SHansol  Suh 
8537f5c9be9SBarry Smith   Input Parameter:
8547f5c9be9SBarry Smith .  tao - the Tao solver context
8557f5c9be9SBarry Smith 
8567f5c9be9SBarry Smith   Output Parameter:
8577f5c9be9SBarry Smith .  reg - the Tao subsolver context
8586285c0a3SHansol  Suh 
8596285c0a3SHansol  Suh   Level: advanced
860a5a2f7acSBarry Smith 
861db781477SPatrick Sanan .seealso: `TAOADMM`
862a5a2f7acSBarry Smith 
8636285c0a3SHansol  Suh @*/
8649371c9d4SSatish Balay PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg) {
8656285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8666285c0a3SHansol  Suh 
8676285c0a3SHansol  Suh   PetscFunctionBegin;
8686285c0a3SHansol  Suh   *reg = am->subsolverZ;
8696285c0a3SHansol  Suh   PetscFunctionReturn(0);
8706285c0a3SHansol  Suh }
8716285c0a3SHansol  Suh 
8726285c0a3SHansol  Suh /*@
8736285c0a3SHansol  Suh   TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM
8746285c0a3SHansol  Suh 
8756285c0a3SHansol  Suh   Collective on Tao
8766285c0a3SHansol  Suh 
8776285c0a3SHansol  Suh   Input Parameters:
8786285c0a3SHansol  Suh + tao - the Tao solver context
8796285c0a3SHansol  Suh - c - RHS vector
8806285c0a3SHansol  Suh 
8816285c0a3SHansol  Suh   Level: advanced
882a5a2f7acSBarry Smith 
883db781477SPatrick Sanan .seealso: `TAOADMM`
884a5a2f7acSBarry Smith 
8856285c0a3SHansol  Suh @*/
8869371c9d4SSatish Balay PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c) {
8876285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8886285c0a3SHansol  Suh 
8896285c0a3SHansol  Suh   PetscFunctionBegin;
8906285c0a3SHansol  Suh   am->constraint = c;
8916285c0a3SHansol  Suh   PetscFunctionReturn(0);
8926285c0a3SHansol  Suh }
8936285c0a3SHansol  Suh 
8946285c0a3SHansol  Suh /*@
8957f5c9be9SBarry Smith   TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty
8966285c0a3SHansol  Suh 
8976285c0a3SHansol  Suh   Collective on Tao
8986285c0a3SHansol  Suh 
8996285c0a3SHansol  Suh   Input Parameters:
9006285c0a3SHansol  Suh +  tao - the Tao solver context
9016285c0a3SHansol  Suh -  mu  - minimum spectral penalty value
9026285c0a3SHansol  Suh 
9036285c0a3SHansol  Suh   Level: advanced
9046285c0a3SHansol  Suh 
905db781477SPatrick Sanan .seealso: `TaoADMMGetSpectralPenalty()`, `TAOADMM`
9066285c0a3SHansol  Suh @*/
9079371c9d4SSatish Balay PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu) {
9086285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
9096285c0a3SHansol  Suh 
9106285c0a3SHansol  Suh   PetscFunctionBegin;
9116285c0a3SHansol  Suh   am->mumin = mu;
9126285c0a3SHansol  Suh   PetscFunctionReturn(0);
9136285c0a3SHansol  Suh }
9146285c0a3SHansol  Suh 
9156285c0a3SHansol  Suh /*@
9166285c0a3SHansol  Suh   TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case
9176285c0a3SHansol  Suh 
9186285c0a3SHansol  Suh   Collective on Tao
9196285c0a3SHansol  Suh 
9206285c0a3SHansol  Suh   Input Parameters:
9216285c0a3SHansol  Suh +  tao - the Tao solver context
9226285c0a3SHansol  Suh -  lambda - L1-norm regularizer coefficient
9236285c0a3SHansol  Suh 
9246285c0a3SHansol  Suh   Level: advanced
9257f5c9be9SBarry Smith 
926db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
9277f5c9be9SBarry Smith 
9286285c0a3SHansol  Suh @*/
9299371c9d4SSatish Balay PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda) {
9306285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
9316285c0a3SHansol  Suh 
9326285c0a3SHansol  Suh   PetscFunctionBegin;
9336285c0a3SHansol  Suh   am->lambda = lambda;
9346285c0a3SHansol  Suh   PetscFunctionReturn(0);
9356285c0a3SHansol  Suh }
9366285c0a3SHansol  Suh 
9376285c0a3SHansol  Suh /*@C
9387f5c9be9SBarry Smith   TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable.
9396285c0a3SHansol  Suh 
9406285c0a3SHansol  Suh   Collective on Tao
9416285c0a3SHansol  Suh 
942a5a2f7acSBarry Smith   Input Parameters:
943a5a2f7acSBarry Smith + tao - the Tao solver context
9446285c0a3SHansol  Suh . J - user-created regularizer constraint Jacobian matrix
9456285c0a3SHansol  Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix
9466285c0a3SHansol  Suh . func - function pointer for the regularizer constraint Jacobian update function
9476285c0a3SHansol  Suh - ctx - user context for the regularizer Hessian
9486285c0a3SHansol  Suh 
9496285c0a3SHansol  Suh   Level: advanced
9507f5c9be9SBarry Smith 
951db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
9527f5c9be9SBarry Smith 
9536285c0a3SHansol  Suh @*/
9549371c9d4SSatish Balay PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx) {
9556285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
9566285c0a3SHansol  Suh 
9576285c0a3SHansol  Suh   PetscFunctionBegin;
9586285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
9596285c0a3SHansol  Suh   if (J) {
9606285c0a3SHansol  Suh     PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
9616285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, J, 2);
9626285c0a3SHansol  Suh   }
9636285c0a3SHansol  Suh   if (Jpre) {
9646285c0a3SHansol  Suh     PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
9656285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Jpre, 3);
9666285c0a3SHansol  Suh   }
9676285c0a3SHansol  Suh   if (ctx) am->misfitjacobianP = ctx;
9686285c0a3SHansol  Suh   if (func) am->ops->misfitjac = func;
9696285c0a3SHansol  Suh 
9706285c0a3SHansol  Suh   if (J) {
9719566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)J));
9729566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JA));
9736285c0a3SHansol  Suh     am->JA = J;
9746285c0a3SHansol  Suh   }
9756285c0a3SHansol  Suh   if (Jpre) {
9769566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Jpre));
9779566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JApre));
9786285c0a3SHansol  Suh     am->JApre = Jpre;
9796285c0a3SHansol  Suh   }
9806285c0a3SHansol  Suh   PetscFunctionReturn(0);
9816285c0a3SHansol  Suh }
9826285c0a3SHansol  Suh 
9836285c0a3SHansol  Suh /*@C
9846285c0a3SHansol  Suh   TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable.
9856285c0a3SHansol  Suh 
9866285c0a3SHansol  Suh   Collective on Tao
9876285c0a3SHansol  Suh 
9886285c0a3SHansol  Suh   Input Parameters:
9896285c0a3SHansol  Suh + tao - the Tao solver context
9906285c0a3SHansol  Suh . J - user-created regularizer constraint Jacobian matrix
9916285c0a3SHansol  Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix
9926285c0a3SHansol  Suh . func - function pointer for the regularizer constraint Jacobian update function
9936285c0a3SHansol  Suh - ctx - user context for the regularizer Hessian
9946285c0a3SHansol  Suh 
9956285c0a3SHansol  Suh   Level: advanced
9967f5c9be9SBarry Smith 
997db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetMisfitConstraintJacobian()`, `TAOADMM`
9987f5c9be9SBarry Smith 
9996285c0a3SHansol  Suh @*/
10009371c9d4SSatish Balay PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx) {
10016285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
10026285c0a3SHansol  Suh 
10036285c0a3SHansol  Suh   PetscFunctionBegin;
10046285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10056285c0a3SHansol  Suh   if (J) {
10066285c0a3SHansol  Suh     PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
10076285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, J, 2);
10086285c0a3SHansol  Suh   }
10096285c0a3SHansol  Suh   if (Jpre) {
10106285c0a3SHansol  Suh     PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
10116285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Jpre, 3);
10126285c0a3SHansol  Suh   }
10136285c0a3SHansol  Suh   if (ctx) am->regjacobianP = ctx;
10146285c0a3SHansol  Suh   if (func) am->ops->regjac = func;
10156285c0a3SHansol  Suh 
10166285c0a3SHansol  Suh   if (J) {
10179566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)J));
10189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JB));
10196285c0a3SHansol  Suh     am->JB = J;
10206285c0a3SHansol  Suh   }
10216285c0a3SHansol  Suh   if (Jpre) {
10229566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Jpre));
10239566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JBpre));
10246285c0a3SHansol  Suh     am->JBpre = Jpre;
10256285c0a3SHansol  Suh   }
10266285c0a3SHansol  Suh   PetscFunctionReturn(0);
10276285c0a3SHansol  Suh }
10286285c0a3SHansol  Suh 
10296285c0a3SHansol  Suh /*@C
10307f5c9be9SBarry Smith    TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function
10317f5c9be9SBarry Smith 
10327f5c9be9SBarry Smith    Collective on tao
10336285c0a3SHansol  Suh 
10346285c0a3SHansol  Suh    Input Parameters:
10356285c0a3SHansol  Suh +    tao - the Tao context
10366285c0a3SHansol  Suh .    func - function pointer for the misfit value and gradient evaluation
10376285c0a3SHansol  Suh -    ctx - user context for the misfit
10386285c0a3SHansol  Suh 
10396285c0a3SHansol  Suh    Level: advanced
10407f5c9be9SBarry Smith 
1041db781477SPatrick Sanan .seealso: `TAOADMM`
10427f5c9be9SBarry Smith 
10436285c0a3SHansol  Suh @*/
10449371c9d4SSatish Balay PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), void *ctx) {
10456285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
10466285c0a3SHansol  Suh 
10476285c0a3SHansol  Suh   PetscFunctionBegin;
10486285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10496285c0a3SHansol  Suh   am->misfitobjgradP     = ctx;
10506285c0a3SHansol  Suh   am->ops->misfitobjgrad = func;
10516285c0a3SHansol  Suh   PetscFunctionReturn(0);
10526285c0a3SHansol  Suh }
10536285c0a3SHansol  Suh 
10546285c0a3SHansol  Suh /*@C
10556285c0a3SHansol  Suh    TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back
10566285c0a3SHansol  Suh    function into the algorithm, to be used for subsolverX.
10576285c0a3SHansol  Suh 
10587f5c9be9SBarry Smith    Collective on tao
10597f5c9be9SBarry Smith 
10606285c0a3SHansol  Suh    Input Parameters:
10616285c0a3SHansol  Suh + tao - the Tao context
10626285c0a3SHansol  Suh . H - user-created matrix for the Hessian of the misfit term
10636285c0a3SHansol  Suh . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term
10646285c0a3SHansol  Suh . func - function pointer for the misfit Hessian evaluation
10656285c0a3SHansol  Suh - ctx - user context for the misfit Hessian
10666285c0a3SHansol  Suh 
10676285c0a3SHansol  Suh    Level: advanced
1068a5a2f7acSBarry Smith 
1069db781477SPatrick Sanan .seealso: `TAOADMM`
1070a5a2f7acSBarry Smith 
10716285c0a3SHansol  Suh @*/
10729371c9d4SSatish Balay PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx) {
10736285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
10746285c0a3SHansol  Suh 
10756285c0a3SHansol  Suh   PetscFunctionBegin;
10766285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10776285c0a3SHansol  Suh   if (H) {
10786285c0a3SHansol  Suh     PetscValidHeaderSpecific(H, MAT_CLASSID, 2);
10796285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, H, 2);
10806285c0a3SHansol  Suh   }
10816285c0a3SHansol  Suh   if (Hpre) {
10826285c0a3SHansol  Suh     PetscValidHeaderSpecific(Hpre, MAT_CLASSID, 3);
10836285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Hpre, 3);
10846285c0a3SHansol  Suh   }
1085ad540459SPierre Jolivet   if (ctx) am->misfithessP = ctx;
1086ad540459SPierre Jolivet   if (func) am->ops->misfithess = func;
10876285c0a3SHansol  Suh   if (H) {
10889566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)H));
10899566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hx));
10906285c0a3SHansol  Suh     am->Hx = H;
10916285c0a3SHansol  Suh   }
10926285c0a3SHansol  Suh   if (Hpre) {
10939566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Hpre));
10949566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hxpre));
10956285c0a3SHansol  Suh     am->Hxpre = Hpre;
10966285c0a3SHansol  Suh   }
10976285c0a3SHansol  Suh   PetscFunctionReturn(0);
10986285c0a3SHansol  Suh }
10996285c0a3SHansol  Suh 
11006285c0a3SHansol  Suh /*@C
11017f5c9be9SBarry Smith    TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function
11027f5c9be9SBarry Smith 
11037f5c9be9SBarry Smith    Collective on tao
11046285c0a3SHansol  Suh 
11056285c0a3SHansol  Suh    Input Parameters:
11066285c0a3SHansol  Suh + tao - the Tao context
11076285c0a3SHansol  Suh . func - function pointer for the regularizer value and gradient evaluation
11086285c0a3SHansol  Suh - ctx - user context for the regularizer
11096285c0a3SHansol  Suh 
11106285c0a3SHansol  Suh    Level: advanced
1111a5a2f7acSBarry Smith 
1112db781477SPatrick Sanan .seealso: `TAOADMM`
1113a5a2f7acSBarry Smith 
11146285c0a3SHansol  Suh @*/
11159371c9d4SSatish Balay PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), void *ctx) {
11166285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
11176285c0a3SHansol  Suh 
11186285c0a3SHansol  Suh   PetscFunctionBegin;
11196285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11206285c0a3SHansol  Suh   am->regobjgradP     = ctx;
11216285c0a3SHansol  Suh   am->ops->regobjgrad = func;
11226285c0a3SHansol  Suh   PetscFunctionReturn(0);
11236285c0a3SHansol  Suh }
11246285c0a3SHansol  Suh 
11256285c0a3SHansol  Suh /*@C
11266285c0a3SHansol  Suh    TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back
11277f5c9be9SBarry Smith    function, to be used for subsolverZ.
11287f5c9be9SBarry Smith 
11297f5c9be9SBarry Smith    Collective on tao
11306285c0a3SHansol  Suh 
11316285c0a3SHansol  Suh    Input Parameters:
11326285c0a3SHansol  Suh + tao - the Tao context
11336285c0a3SHansol  Suh . H - user-created matrix for the Hessian of the regularization term
11346285c0a3SHansol  Suh . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term
11356285c0a3SHansol  Suh . func - function pointer for the regularizer Hessian evaluation
11366285c0a3SHansol  Suh - ctx - user context for the regularizer Hessian
11376285c0a3SHansol  Suh 
11386285c0a3SHansol  Suh    Level: advanced
1139a5a2f7acSBarry Smith 
1140db781477SPatrick Sanan .seealso: `TAOADMM`
1141a5a2f7acSBarry Smith 
11426285c0a3SHansol  Suh @*/
11439371c9d4SSatish Balay PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx) {
11446285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
11456285c0a3SHansol  Suh 
11466285c0a3SHansol  Suh   PetscFunctionBegin;
11476285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11486285c0a3SHansol  Suh   if (H) {
11496285c0a3SHansol  Suh     PetscValidHeaderSpecific(H, MAT_CLASSID, 2);
11506285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, H, 2);
11516285c0a3SHansol  Suh   }
11526285c0a3SHansol  Suh   if (Hpre) {
11536285c0a3SHansol  Suh     PetscValidHeaderSpecific(Hpre, MAT_CLASSID, 3);
11546285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Hpre, 3);
11556285c0a3SHansol  Suh   }
1156ad540459SPierre Jolivet   if (ctx) am->reghessP = ctx;
1157ad540459SPierre Jolivet   if (func) am->ops->reghess = func;
11586285c0a3SHansol  Suh   if (H) {
11599566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)H));
11609566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hz));
11616285c0a3SHansol  Suh     am->Hz = H;
11626285c0a3SHansol  Suh   }
11636285c0a3SHansol  Suh   if (Hpre) {
11649566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Hpre));
11659566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hzpre));
11666285c0a3SHansol  Suh     am->Hzpre = Hpre;
11676285c0a3SHansol  Suh   }
11686285c0a3SHansol  Suh   PetscFunctionReturn(0);
11696285c0a3SHansol  Suh }
11706285c0a3SHansol  Suh 
11716285c0a3SHansol  Suh /*@
11727f5c9be9SBarry Smith    TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver.
11736285c0a3SHansol  Suh 
11747f5c9be9SBarry Smith    Collective on tao
11757f5c9be9SBarry Smith 
11767f5c9be9SBarry Smith    Input Parameter:
11777f5c9be9SBarry Smith . tao - the Tao context
11787f5c9be9SBarry Smith 
11797f5c9be9SBarry Smith    Output Parameter:
11807f5c9be9SBarry Smith . admm_tao - the parent Tao context
11816285c0a3SHansol  Suh 
11826285c0a3SHansol  Suh    Level: advanced
1183a5a2f7acSBarry Smith 
1184db781477SPatrick Sanan .seealso: `TAOADMM`
1185a5a2f7acSBarry Smith 
11866285c0a3SHansol  Suh @*/
11879371c9d4SSatish Balay PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao) {
11886285c0a3SHansol  Suh   PetscFunctionBegin;
11896285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11909566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)tao, "TaoGetADMMParentTao_ADMM", (PetscObject *)admm_tao));
11916285c0a3SHansol  Suh   PetscFunctionReturn(0);
11926285c0a3SHansol  Suh }
11936285c0a3SHansol  Suh 
11946285c0a3SHansol  Suh /*@
11957f5c9be9SBarry Smith   TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state
11966285c0a3SHansol  Suh 
11976285c0a3SHansol  Suh   Not Collective
11986285c0a3SHansol  Suh 
11996285c0a3SHansol  Suh   Input Parameter:
12007f5c9be9SBarry Smith . tao - the Tao context
12017f5c9be9SBarry Smith 
12027f5c9be9SBarry Smith   Output Parameter:
12037f5c9be9SBarry Smith . Y - the current solution
12046285c0a3SHansol  Suh 
12056285c0a3SHansol  Suh   Level: intermediate
1206a5a2f7acSBarry Smith 
1207db781477SPatrick Sanan .seealso: `TAOADMM`
1208a5a2f7acSBarry Smith 
12096285c0a3SHansol  Suh @*/
12109371c9d4SSatish Balay PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y) {
12116285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
12126285c0a3SHansol  Suh 
12136285c0a3SHansol  Suh   PetscFunctionBegin;
12146285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12156285c0a3SHansol  Suh   *Y = am->y;
12166285c0a3SHansol  Suh   PetscFunctionReturn(0);
12176285c0a3SHansol  Suh }
12186285c0a3SHansol  Suh 
12196285c0a3SHansol  Suh /*@
12206285c0a3SHansol  Suh   TaoADMMSetRegularizerType - Set regularizer type for ADMM routine
12216285c0a3SHansol  Suh 
12226285c0a3SHansol  Suh   Not Collective
12236285c0a3SHansol  Suh 
1224d8d19677SJose E. Roman   Input Parameters:
12256285c0a3SHansol  Suh + tao  - the Tao context
12267f5c9be9SBarry Smith - type - regularizer type
12277f5c9be9SBarry Smith 
12287f5c9be9SBarry Smith   Options Database:
1229147403d9SBarry Smith .  -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer
12306285c0a3SHansol  Suh 
12316285c0a3SHansol  Suh   Level: intermediate
12326285c0a3SHansol  Suh 
1233db781477SPatrick Sanan .seealso: `TaoADMMGetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
12346285c0a3SHansol  Suh @*/
12359371c9d4SSatish Balay PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type) {
12366285c0a3SHansol  Suh   PetscFunctionBegin;
12376285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12386285c0a3SHansol  Suh   PetscValidLogicalCollectiveEnum(tao, type, 2);
1239cac4c232SBarry Smith   PetscTryMethod(tao, "TaoADMMSetRegularizerType_C", (Tao, TaoADMMRegularizerType), (tao, type));
12406285c0a3SHansol  Suh   PetscFunctionReturn(0);
12416285c0a3SHansol  Suh }
12426285c0a3SHansol  Suh 
12436285c0a3SHansol  Suh /*@
12446285c0a3SHansol  Suh    TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM
12456285c0a3SHansol  Suh 
12466285c0a3SHansol  Suh    Not Collective
12476285c0a3SHansol  Suh 
12486285c0a3SHansol  Suh    Input Parameter:
12496285c0a3SHansol  Suh .  tao - the Tao context
12506285c0a3SHansol  Suh 
12516285c0a3SHansol  Suh    Output Parameter:
12526285c0a3SHansol  Suh .  type - the type of regularizer
12536285c0a3SHansol  Suh 
12546285c0a3SHansol  Suh    Level: intermediate
12556285c0a3SHansol  Suh 
1256db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
12576285c0a3SHansol  Suh @*/
12589371c9d4SSatish Balay PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type) {
12596285c0a3SHansol  Suh   PetscFunctionBegin;
12606285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1261cac4c232SBarry Smith   PetscUseMethod(tao, "TaoADMMGetRegularizerType_C", (Tao, TaoADMMRegularizerType *), (tao, type));
12626285c0a3SHansol  Suh   PetscFunctionReturn(0);
12636285c0a3SHansol  Suh }
12646285c0a3SHansol  Suh 
12656285c0a3SHansol  Suh /*@
12666285c0a3SHansol  Suh   TaoADMMSetUpdateType - Set update routine for ADMM routine
12676285c0a3SHansol  Suh 
12686285c0a3SHansol  Suh   Not Collective
12696285c0a3SHansol  Suh 
1270d8d19677SJose E. Roman   Input Parameters:
12716285c0a3SHansol  Suh + tao  - the Tao context
12726285c0a3SHansol  Suh - type - spectral parameter update type
12736285c0a3SHansol  Suh 
12746285c0a3SHansol  Suh   Level: intermediate
12756285c0a3SHansol  Suh 
1276db781477SPatrick Sanan .seealso: `TaoADMMGetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
12776285c0a3SHansol  Suh @*/
12789371c9d4SSatish Balay PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type) {
12796285c0a3SHansol  Suh   PetscFunctionBegin;
12806285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12816285c0a3SHansol  Suh   PetscValidLogicalCollectiveEnum(tao, type, 2);
1282cac4c232SBarry Smith   PetscTryMethod(tao, "TaoADMMSetUpdateType_C", (Tao, TaoADMMUpdateType), (tao, type));
12836285c0a3SHansol  Suh   PetscFunctionReturn(0);
12846285c0a3SHansol  Suh }
12856285c0a3SHansol  Suh 
12866285c0a3SHansol  Suh /*@
12876285c0a3SHansol  Suh    TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM
12886285c0a3SHansol  Suh 
12896285c0a3SHansol  Suh    Not Collective
12906285c0a3SHansol  Suh 
12916285c0a3SHansol  Suh    Input Parameter:
1292f0fc11ceSJed Brown .  tao - the Tao context
12936285c0a3SHansol  Suh 
12946285c0a3SHansol  Suh    Output Parameter:
12956285c0a3SHansol  Suh .  type - the type of spectral penalty update routine
12966285c0a3SHansol  Suh 
12976285c0a3SHansol  Suh    Level: intermediate
12986285c0a3SHansol  Suh 
1299db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
13006285c0a3SHansol  Suh @*/
13019371c9d4SSatish Balay PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type) {
13026285c0a3SHansol  Suh   PetscFunctionBegin;
13036285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1304cac4c232SBarry Smith   PetscUseMethod(tao, "TaoADMMGetUpdateType_C", (Tao, TaoADMMUpdateType *), (tao, type));
13056285c0a3SHansol  Suh   PetscFunctionReturn(0);
13066285c0a3SHansol  Suh }
1307