xref: /petsc/src/tao/constrained/impls/admm/admm.c (revision b4623fece54afb61fe62f4c6eafc62dfcda9fd50)
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 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoADMMToleranceUpdate(Tao tao)
28d71ae5a4SJacob Faibussowitsch {
296285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
306285c0a3SHansol  Suh   PetscReal Axnorm, Bznorm, ATynorm, temp;
316285c0a3SHansol  Suh   Vec       tempJR, tempL;
326285c0a3SHansol  Suh   Tao       mis;
336285c0a3SHansol  Suh 
346285c0a3SHansol  Suh   PetscFunctionBegin;
356285c0a3SHansol  Suh   mis    = am->subsolverX;
366285c0a3SHansol  Suh   tempJR = am->workJacobianRight;
376285c0a3SHansol  Suh   tempL  = am->workLeft;
386285c0a3SHansol  Suh   /* ATy */
399566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianEquality(mis, am->y, mis->jacobian_equality, mis->jacobian_equality_pre));
409566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mis->jacobian_equality, am->y, tempJR));
419566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR, NORM_2, &ATynorm));
426285c0a3SHansol  Suh   /* dualres = mu * ||AT(Bz-Bzold)||_2 */
439566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR, -1., am->Bzold, am->Bz));
449566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mis->jacobian_equality, tempJR, tempL));
459566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempL, NORM_2, &am->dualres));
466285c0a3SHansol  Suh   am->dualres *= am->mu;
476285c0a3SHansol  Suh 
486285c0a3SHansol  Suh   /* ||Ax||_2, ||Bz||_2 */
499566063dSJacob Faibussowitsch   PetscCall(VecNorm(am->Ax, NORM_2, &Axnorm));
509566063dSJacob Faibussowitsch   PetscCall(VecNorm(am->Bz, NORM_2, &Bznorm));
516285c0a3SHansol  Suh 
526285c0a3SHansol  Suh   /* Set catol to be catol_admm *  max{||Ax||,||Bz||,||c||} *
536285c0a3SHansol  Suh    * Set gatol to be gatol_admm *  ||A^Ty|| *
546285c0a3SHansol  Suh    * while cnorm is ||r_k||_2, and gnorm is ||d_k||_2 */
556285c0a3SHansol  Suh   temp = am->catol_admm * PetscMax(Axnorm, (!am->const_norm) ? Bznorm : PetscMax(Bznorm, am->const_norm));
569566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintTolerances(tao, temp, PETSC_DEFAULT));
579566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao, am->gatol_admm * ATynorm, PETSC_DEFAULT, PETSC_DEFAULT));
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
596285c0a3SHansol  Suh }
606285c0a3SHansol  Suh 
616285c0a3SHansol  Suh /* Penaly Update for Adaptive ADMM. */
62d71ae5a4SJacob Faibussowitsch static PetscErrorCode AdaptiveADMMPenaltyUpdate(Tao tao)
63d71ae5a4SJacob Faibussowitsch {
646285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
656285c0a3SHansol  Suh   PetscReal ydiff_norm, yhatdiff_norm, Axdiff_norm, Bzdiff_norm, Axyhat, Bzy, a_sd, a_mg, a_k, b_sd, b_mg, b_k;
666285c0a3SHansol  Suh   PetscBool hflag, gflag;
676285c0a3SHansol  Suh   Vec       tempJR, tempJR2;
686285c0a3SHansol  Suh 
696285c0a3SHansol  Suh   PetscFunctionBegin;
706285c0a3SHansol  Suh   tempJR  = am->workJacobianRight;
716285c0a3SHansol  Suh   tempJR2 = am->workJacobianRight2;
726285c0a3SHansol  Suh   hflag   = PETSC_FALSE;
736285c0a3SHansol  Suh   gflag   = PETSC_FALSE;
746285c0a3SHansol  Suh   a_k     = -1;
756285c0a3SHansol  Suh   b_k     = -1;
766285c0a3SHansol  Suh 
779566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR, -1., am->Axold, am->Ax));
789566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR2, -1., am->yhatold, am->yhat));
799566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR, NORM_2, &Axdiff_norm));
809566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR2, NORM_2, &yhatdiff_norm));
819566063dSJacob Faibussowitsch   PetscCall(VecDot(tempJR, tempJR2, &Axyhat));
826285c0a3SHansol  Suh 
839566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR, -1., am->Bz0, am->Bz));
849566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tempJR2, -1., am->y, am->y0));
859566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR, NORM_2, &Bzdiff_norm));
869566063dSJacob Faibussowitsch   PetscCall(VecNorm(tempJR2, NORM_2, &ydiff_norm));
879566063dSJacob Faibussowitsch   PetscCall(VecDot(tempJR, tempJR2, &Bzy));
886285c0a3SHansol  Suh 
896285c0a3SHansol  Suh   if (Axyhat > am->orthval * Axdiff_norm * yhatdiff_norm + am->mueps) {
906285c0a3SHansol  Suh     hflag = PETSC_TRUE;
916285c0a3SHansol  Suh     a_sd  = PetscSqr(yhatdiff_norm) / Axyhat; /* alphaSD */
926285c0a3SHansol  Suh     a_mg  = Axyhat / PetscSqr(Axdiff_norm);   /* alphaMG */
936285c0a3SHansol  Suh     a_k   = (a_mg / a_sd) > 0.5 ? a_mg : a_sd - 0.5 * a_mg;
946285c0a3SHansol  Suh   }
956285c0a3SHansol  Suh   if (Bzy > am->orthval * Bzdiff_norm * ydiff_norm + am->mueps) {
966285c0a3SHansol  Suh     gflag = PETSC_TRUE;
976285c0a3SHansol  Suh     b_sd  = PetscSqr(ydiff_norm) / Bzy;  /* betaSD */
986285c0a3SHansol  Suh     b_mg  = Bzy / PetscSqr(Bzdiff_norm); /* betaMG */
996285c0a3SHansol  Suh     b_k   = (b_mg / b_sd) > 0.5 ? b_mg : b_sd - 0.5 * b_mg;
1006285c0a3SHansol  Suh   }
1016285c0a3SHansol  Suh   am->muold = am->mu;
1026285c0a3SHansol  Suh   if (gflag && hflag) {
1036285c0a3SHansol  Suh     am->mu = PetscSqrtReal(a_k * b_k);
1046285c0a3SHansol  Suh   } else if (hflag) {
1056285c0a3SHansol  Suh     am->mu = a_k;
1066285c0a3SHansol  Suh   } else if (gflag) {
1076285c0a3SHansol  Suh     am->mu = b_k;
1086285c0a3SHansol  Suh   }
109ad540459SPierre Jolivet   if (am->mu > am->muold) am->mu = am->muold;
110ad540459SPierre Jolivet   if (am->mu < am->mumin) am->mu = am->mumin;
1113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1126285c0a3SHansol  Suh }
1136285c0a3SHansol  Suh 
114d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoADMMSetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType type)
115d71ae5a4SJacob Faibussowitsch {
1166285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1176285c0a3SHansol  Suh 
1186285c0a3SHansol  Suh   PetscFunctionBegin;
1196285c0a3SHansol  Suh   am->regswitch = type;
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1216285c0a3SHansol  Suh }
1226285c0a3SHansol  Suh 
123d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoADMMGetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType *type)
124d71ae5a4SJacob Faibussowitsch {
1256285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1266285c0a3SHansol  Suh 
1276285c0a3SHansol  Suh   PetscFunctionBegin;
1286285c0a3SHansol  Suh   *type = am->regswitch;
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1306285c0a3SHansol  Suh }
1316285c0a3SHansol  Suh 
132d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoADMMSetUpdateType_ADMM(Tao tao, TaoADMMUpdateType type)
133d71ae5a4SJacob Faibussowitsch {
1346285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1356285c0a3SHansol  Suh 
1366285c0a3SHansol  Suh   PetscFunctionBegin;
1376285c0a3SHansol  Suh   am->update = type;
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1396285c0a3SHansol  Suh }
1406285c0a3SHansol  Suh 
141d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoADMMGetUpdateType_ADMM(Tao tao, TaoADMMUpdateType *type)
142d71ae5a4SJacob Faibussowitsch {
1436285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1446285c0a3SHansol  Suh 
1456285c0a3SHansol  Suh   PetscFunctionBegin;
1466285c0a3SHansol  Suh   *type = am->update;
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1486285c0a3SHansol  Suh }
1496285c0a3SHansol  Suh 
1506285c0a3SHansol  Suh /* This routine updates Jacobians with new x,z vectors,
1516285c0a3SHansol  Suh  * and then updates Ax and Bz vectors, then computes updated residual vector*/
152d71ae5a4SJacob Faibussowitsch static PetscErrorCode ADMMUpdateConstraintResidualVector(Tao tao, Vec x, Vec z, Vec Ax, Vec Bz, Vec residual)
153d71ae5a4SJacob Faibussowitsch {
1546285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
1556285c0a3SHansol  Suh   Tao       mis, reg;
1566285c0a3SHansol  Suh 
1576285c0a3SHansol  Suh   PetscFunctionBegin;
1586285c0a3SHansol  Suh   mis = am->subsolverX;
1596285c0a3SHansol  Suh   reg = am->subsolverZ;
1609566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianEquality(mis, x, mis->jacobian_equality, mis->jacobian_equality_pre));
1619566063dSJacob Faibussowitsch   PetscCall(MatMult(mis->jacobian_equality, x, Ax));
1629566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianEquality(reg, z, reg->jacobian_equality, reg->jacobian_equality_pre));
1639566063dSJacob Faibussowitsch   PetscCall(MatMult(reg->jacobian_equality, z, Bz));
1646285c0a3SHansol  Suh 
1659566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(residual, 1., Bz, Ax));
16648a46eb9SPierre Jolivet   if (am->constraint != NULL) PetscCall(VecAXPY(residual, -1., am->constraint));
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1686285c0a3SHansol  Suh }
1696285c0a3SHansol  Suh 
1706285c0a3SHansol  Suh /* Updates Augmented Lagrangians to given routines *
1716285c0a3SHansol  Suh  * For subsolverX, routine needs to be ComputeObjectiveAndGraidnet
1726285c0a3SHansol  Suh  * Separate Objective and Gradient routines are not supported.  */
173d71ae5a4SJacob Faibussowitsch static PetscErrorCode SubObjGradUpdate(Tao tao, Vec x, PetscReal *f, Vec g, void *ptr)
174d71ae5a4SJacob Faibussowitsch {
1756285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
1766285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
1776285c0a3SHansol  Suh   PetscReal temp, temp2;
1786285c0a3SHansol  Suh   Vec       tempJR;
1796285c0a3SHansol  Suh 
1806285c0a3SHansol  Suh   PetscFunctionBegin;
1816285c0a3SHansol  Suh   tempJR = am->workJacobianRight;
1829566063dSJacob Faibussowitsch   PetscCall(ADMMUpdateConstraintResidualVector(parent, x, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
1839566063dSJacob Faibussowitsch   PetscCall((*am->ops->misfitobjgrad)(am->subsolverX, x, f, g, am->misfitobjgradP));
1846285c0a3SHansol  Suh 
1856285c0a3SHansol  Suh   am->last_misfit_val = *f;
1866285c0a3SHansol  Suh   /* Objective  Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */
1879566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->y, &temp));
1889566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->residual, &temp2));
1896285c0a3SHansol  Suh   *f += temp + (am->mu / 2) * temp2;
1906285c0a3SHansol  Suh 
1916285c0a3SHansol  Suh   /* Gradient. Add + mu*AT(Ax+Bz-c) + yTA*/
1929566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian_equality, am->residual, tempJR));
1939566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, am->mu, tempJR));
1949566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian_equality, am->y, tempJR));
1959566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, 1., tempJR));
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1976285c0a3SHansol  Suh }
1986285c0a3SHansol  Suh 
1996285c0a3SHansol  Suh /* Updates Augmented Lagrangians to given routines
2006285c0a3SHansol  Suh  * For subsolverZ, routine needs to be ComputeObjectiveAndGraidnet
2016285c0a3SHansol  Suh  * Separate Objective and Gradient routines are not supported.  */
202d71ae5a4SJacob Faibussowitsch static PetscErrorCode RegObjGradUpdate(Tao tao, Vec z, PetscReal *f, Vec g, void *ptr)
203d71ae5a4SJacob Faibussowitsch {
2046285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
2056285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
2066285c0a3SHansol  Suh   PetscReal temp, temp2;
2076285c0a3SHansol  Suh   Vec       tempJR;
2086285c0a3SHansol  Suh 
2096285c0a3SHansol  Suh   PetscFunctionBegin;
2106285c0a3SHansol  Suh   tempJR = am->workJacobianRight;
2119566063dSJacob Faibussowitsch   PetscCall(ADMMUpdateConstraintResidualVector(parent, am->subsolverX->solution, z, am->Ax, am->Bz, am->residual));
2129566063dSJacob Faibussowitsch   PetscCall((*am->ops->regobjgrad)(am->subsolverZ, z, f, g, am->regobjgradP));
2136285c0a3SHansol  Suh   am->last_reg_val = *f;
2146285c0a3SHansol  Suh   /* Objective  Add  + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */
2159566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->y, &temp));
2169566063dSJacob Faibussowitsch   PetscCall(VecTDot(am->residual, am->residual, &temp2));
2176285c0a3SHansol  Suh   *f += temp + (am->mu / 2) * temp2;
2186285c0a3SHansol  Suh 
2196285c0a3SHansol  Suh   /* Gradient. Add + mu*BT(Ax+Bz-c) + yTB*/
2209566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality, am->residual, tempJR));
2219566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, am->mu, tempJR));
2229566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality, am->y, tempJR));
2239566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, 1., tempJR));
2243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2256285c0a3SHansol  Suh }
2266285c0a3SHansol  Suh 
2276285c0a3SHansol  Suh /* Computes epsilon padded L1 norm lambda*sum(sqrt(x^2+eps^2)-eps */
228d71ae5a4SJacob Faibussowitsch static PetscErrorCode ADMML1EpsilonNorm(Tao tao, Vec x, PetscReal eps, PetscReal *norm)
229d71ae5a4SJacob Faibussowitsch {
2306285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
2316285c0a3SHansol  Suh   PetscInt  N;
2326285c0a3SHansol  Suh 
2336285c0a3SHansol  Suh   PetscFunctionBegin;
2349566063dSJacob Faibussowitsch   PetscCall(VecGetSize(am->workLeft, &N));
2359566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(am->workLeft, x, x));
2369566063dSJacob Faibussowitsch   PetscCall(VecShift(am->workLeft, am->l1epsilon * am->l1epsilon));
2379566063dSJacob Faibussowitsch   PetscCall(VecSqrtAbs(am->workLeft));
2389566063dSJacob Faibussowitsch   PetscCall(VecSum(am->workLeft, norm));
2396285c0a3SHansol  Suh   *norm += N * am->l1epsilon;
2406285c0a3SHansol  Suh   *norm *= am->lambda;
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2426285c0a3SHansol  Suh }
2436285c0a3SHansol  Suh 
244d71ae5a4SJacob Faibussowitsch static PetscErrorCode ADMMInternalHessianUpdate(Mat H, Mat Constraint, PetscBool Identity, void *ptr)
245d71ae5a4SJacob Faibussowitsch {
2466285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)ptr;
2476285c0a3SHansol  Suh 
2486285c0a3SHansol  Suh   PetscFunctionBegin;
2496285c0a3SHansol  Suh   switch (am->update) {
250d71ae5a4SJacob Faibussowitsch   case (TAO_ADMM_UPDATE_BASIC):
251d71ae5a4SJacob Faibussowitsch     break;
2526285c0a3SHansol  Suh   case (TAO_ADMM_UPDATE_ADAPTIVE):
2536285c0a3SHansol  Suh   case (TAO_ADMM_UPDATE_ADAPTIVE_RELAXED):
2546285c0a3SHansol  Suh     if (H && (am->muold != am->mu)) {
2556285c0a3SHansol  Suh       if (!Identity) {
2569566063dSJacob Faibussowitsch         PetscCall(MatAXPY(H, am->mu - am->muold, Constraint, DIFFERENT_NONZERO_PATTERN));
2576285c0a3SHansol  Suh       } else {
2589566063dSJacob Faibussowitsch         PetscCall(MatShift(H, am->mu - am->muold));
2596285c0a3SHansol  Suh       }
2606285c0a3SHansol  Suh     }
2616285c0a3SHansol  Suh     break;
2626285c0a3SHansol  Suh   }
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2646285c0a3SHansol  Suh }
2656285c0a3SHansol  Suh 
2666285c0a3SHansol  Suh /* Updates Hessian - adds second derivative of augmented Lagrangian
2676285c0a3SHansol  Suh  * H \gets H + \rho*ATA
26869d47153SPierre Jolivet   Here, \rho does not change in TAO_ADMM_UPDATE_BASIC - thus no-op
26969d47153SPierre Jolivet   For ADAPTAIVE,ADAPTIVE_RELAXED,
27069d47153SPierre Jolivet   H \gets H + (\rho-\rhoold)*ATA
27169d47153SPierre Jolivet   Here, we assume that A is linear constraint i.e., does not change.
27269d47153SPierre Jolivet   Thus, for both ADAPTIVE, and RELAXED, ATA matrix is pre-set (except for A=I (null case)) see TaoSetUp_ADMM */
273d71ae5a4SJacob Faibussowitsch static PetscErrorCode SubHessianUpdate(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
274d71ae5a4SJacob Faibussowitsch {
2756285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
2766285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
2776285c0a3SHansol  Suh 
2786285c0a3SHansol  Suh   PetscFunctionBegin;
2796285c0a3SHansol  Suh   if (am->Hxchange) {
2806285c0a3SHansol  Suh     /* Case where Hessian gets updated with respect to x vector input. */
2819566063dSJacob Faibussowitsch     PetscCall((*am->ops->misfithess)(am->subsolverX, x, H, Hpre, am->misfithessP));
2829566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian, am->ATA, am->xJI, am));
2836285c0a3SHansol  Suh   } else if (am->Hxbool) {
2846285c0a3SHansol  Suh     /* Hessian doesn't get updated. H(x) = c */
2856285c0a3SHansol  Suh     /* Update Lagrangian only only per TAO call */
2869566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian, am->ATA, am->xJI, am));
2876285c0a3SHansol  Suh     am->Hxbool = PETSC_FALSE;
2886285c0a3SHansol  Suh   }
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2906285c0a3SHansol  Suh }
2916285c0a3SHansol  Suh 
2926285c0a3SHansol  Suh /* Same as SubHessianUpdate, except for B matrix instead of A matrix */
293d71ae5a4SJacob Faibussowitsch static PetscErrorCode RegHessianUpdate(Tao tao, Vec z, Mat H, Mat Hpre, void *ptr)
294d71ae5a4SJacob Faibussowitsch {
2956285c0a3SHansol  Suh   Tao       parent = (Tao)ptr;
2966285c0a3SHansol  Suh   TAO_ADMM *am     = (TAO_ADMM *)parent->data;
2976285c0a3SHansol  Suh 
2986285c0a3SHansol  Suh   PetscFunctionBegin;
2996285c0a3SHansol  Suh 
3006285c0a3SHansol  Suh   if (am->Hzchange) {
3016285c0a3SHansol  Suh     /* Case where Hessian gets updated with respect to x vector input. */
3029566063dSJacob Faibussowitsch     PetscCall((*am->ops->reghess)(am->subsolverZ, z, H, Hpre, am->reghessP));
3039566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian, am->BTB, am->zJI, am));
3046285c0a3SHansol  Suh   } else if (am->Hzbool) {
3056285c0a3SHansol  Suh     /* Hessian doesn't get updated. H(x) = c */
3066285c0a3SHansol  Suh     /* Update Lagrangian only only per TAO call */
3079566063dSJacob Faibussowitsch     PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian, am->BTB, am->zJI, am));
3086285c0a3SHansol  Suh     am->Hzbool = PETSC_FALSE;
3096285c0a3SHansol  Suh   }
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3116285c0a3SHansol  Suh }
3126285c0a3SHansol  Suh 
3136285c0a3SHansol  Suh /* Shell Matrix routine for A matrix.
3146285c0a3SHansol  Suh  * This gets used when user puts NULL for
3156285c0a3SHansol  Suh  * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...)
3166285c0a3SHansol  Suh  * Essentially sets A=I*/
317d71ae5a4SJacob Faibussowitsch static PetscErrorCode JacobianIdentity(Mat mat, Vec in, Vec out)
318d71ae5a4SJacob Faibussowitsch {
3196285c0a3SHansol  Suh   PetscFunctionBegin;
3209566063dSJacob Faibussowitsch   PetscCall(VecCopy(in, out));
3213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3226285c0a3SHansol  Suh }
3236285c0a3SHansol  Suh 
3246285c0a3SHansol  Suh /* Shell Matrix routine for B matrix.
3256285c0a3SHansol  Suh  * This gets used when user puts NULL for
3266285c0a3SHansol  Suh  * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...)
3276285c0a3SHansol  Suh  * Sets B=-I */
328d71ae5a4SJacob Faibussowitsch static PetscErrorCode JacobianIdentityB(Mat mat, Vec in, Vec out)
329d71ae5a4SJacob Faibussowitsch {
3306285c0a3SHansol  Suh   PetscFunctionBegin;
3319566063dSJacob Faibussowitsch   PetscCall(VecCopy(in, out));
3329566063dSJacob Faibussowitsch   PetscCall(VecScale(out, -1.));
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3346285c0a3SHansol  Suh }
3356285c0a3SHansol  Suh 
3366285c0a3SHansol  Suh /* Solve f(x) + g(z) s.t. Ax + Bz = c */
337d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_ADMM(Tao tao)
338d71ae5a4SJacob Faibussowitsch {
3396285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
3406285c0a3SHansol  Suh   PetscInt  N;
3416285c0a3SHansol  Suh   PetscReal reg_func;
3426285c0a3SHansol  Suh   PetscBool is_reg_shell;
3436285c0a3SHansol  Suh   Vec       tempL;
3446285c0a3SHansol  Suh 
3456285c0a3SHansol  Suh   PetscFunctionBegin;
3466285c0a3SHansol  Suh   if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) {
3473c859ba3SBarry Smith     PetscCheck(am->subsolverX->ops->computejacobianequality, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoADMMSetMisfitConstraintJacobian() first");
3483c859ba3SBarry Smith     PetscCheck(am->subsolverZ->ops->computejacobianequality, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoADMMSetRegularizerConstraintJacobian() first");
34948a46eb9SPierre Jolivet     if (am->constraint != NULL) PetscCall(VecNorm(am->constraint, NORM_2, &am->const_norm));
3506285c0a3SHansol  Suh   }
3516285c0a3SHansol  Suh   tempL = am->workLeft;
3529566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tempL, &N));
3536285c0a3SHansol  Suh 
35448a46eb9SPierre Jolivet   if (am->Hx && am->ops->misfithess) PetscCall(TaoSetHessian(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao));
3556285c0a3SHansol  Suh 
3566285c0a3SHansol  Suh   if (!am->zJI) {
3576285c0a3SHansol  Suh     /* Currently, B is assumed to be a linear system, i.e., not getting updated*/
3589566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(am->JB, am->JB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(am->BTB)));
3596285c0a3SHansol  Suh   }
3606285c0a3SHansol  Suh   if (!am->xJI) {
3616285c0a3SHansol  Suh     /* Currently, A is assumed to be a linear system, i.e., not getting updated*/
3629566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(am->ATA)));
3636285c0a3SHansol  Suh   }
3646285c0a3SHansol  Suh 
3656285c0a3SHansol  Suh   is_reg_shell = PETSC_FALSE;
3666285c0a3SHansol  Suh 
3679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell));
3686285c0a3SHansol  Suh 
3696285c0a3SHansol  Suh   if (!is_reg_shell) {
3706285c0a3SHansol  Suh     switch (am->regswitch) {
371d71ae5a4SJacob Faibussowitsch     case (TAO_ADMM_REGULARIZER_USER):
372d71ae5a4SJacob Faibussowitsch       break;
3736285c0a3SHansol  Suh     case (TAO_ADMM_REGULARIZER_SOFT_THRESH):
3746285c0a3SHansol  Suh       /* Soft Threshold. */
3756285c0a3SHansol  Suh       break;
3766285c0a3SHansol  Suh     }
3771baa6e33SBarry Smith     if (am->ops->regobjgrad) PetscCall(TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao));
37848a46eb9SPierre Jolivet     if (am->Hz && am->ops->reghess) PetscCall(TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao));
3796285c0a3SHansol  Suh   }
3806285c0a3SHansol  Suh 
3816285c0a3SHansol  Suh   switch (am->update) {
3826285c0a3SHansol  Suh   case TAO_ADMM_UPDATE_BASIC:
3836285c0a3SHansol  Suh     if (am->subsolverX->hessian) {
3846285c0a3SHansol  Suh       /* In basic case, Hessian does not get updated w.r.t. to spectral penalty
3856285c0a3SHansol  Suh        * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/
3866285c0a3SHansol  Suh       if (!am->xJI) {
3879566063dSJacob Faibussowitsch         PetscCall(MatAXPY(am->subsolverX->hessian, am->mu, am->ATA, DIFFERENT_NONZERO_PATTERN));
3886285c0a3SHansol  Suh       } else {
3899566063dSJacob Faibussowitsch         PetscCall(MatShift(am->subsolverX->hessian, am->mu));
3906285c0a3SHansol  Suh       }
3916285c0a3SHansol  Suh     }
3926285c0a3SHansol  Suh     if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) {
3936285c0a3SHansol  Suh       if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) {
3949566063dSJacob Faibussowitsch         PetscCall(MatAXPY(am->subsolverZ->hessian, am->mu, am->BTB, DIFFERENT_NONZERO_PATTERN));
3956285c0a3SHansol  Suh       } else {
3969566063dSJacob Faibussowitsch         PetscCall(MatShift(am->subsolverZ->hessian, am->mu));
3976285c0a3SHansol  Suh       }
3986285c0a3SHansol  Suh     }
3996285c0a3SHansol  Suh     break;
4006285c0a3SHansol  Suh   case TAO_ADMM_UPDATE_ADAPTIVE:
401d71ae5a4SJacob Faibussowitsch   case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
402d71ae5a4SJacob Faibussowitsch     break;
4036285c0a3SHansol  Suh   }
4046285c0a3SHansol  Suh 
4059566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
4066285c0a3SHansol  Suh   tao->reason = TAO_CONTINUE_ITERATING;
4076285c0a3SHansol  Suh 
4086285c0a3SHansol  Suh   while (tao->reason == TAO_CONTINUE_ITERATING) {
409dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
4109566063dSJacob Faibussowitsch     PetscCall(VecCopy(am->Bz, am->Bzold));
4116285c0a3SHansol  Suh 
4126285c0a3SHansol  Suh     /* x update */
4139566063dSJacob Faibussowitsch     PetscCall(TaoSolve(am->subsolverX));
4149566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre));
4159566063dSJacob Faibussowitsch     PetscCall(MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution, am->Ax));
4166285c0a3SHansol  Suh 
4176285c0a3SHansol  Suh     am->Hxbool = PETSC_TRUE;
4186285c0a3SHansol  Suh 
4196285c0a3SHansol  Suh     /* z update */
4206285c0a3SHansol  Suh     switch (am->regswitch) {
421d71ae5a4SJacob Faibussowitsch     case TAO_ADMM_REGULARIZER_USER:
422d71ae5a4SJacob Faibussowitsch       PetscCall(TaoSolve(am->subsolverZ));
423d71ae5a4SJacob Faibussowitsch       break;
4246285c0a3SHansol  Suh     case TAO_ADMM_REGULARIZER_SOFT_THRESH:
4256285c0a3SHansol  Suh       /* L1 assumes A,B jacobians are identity nxn matrix */
4269566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(am->workJacobianRight, 1 / am->mu, am->y, am->Ax));
4279566063dSJacob Faibussowitsch       PetscCall(TaoSoftThreshold(am->workJacobianRight, -am->lambda / am->mu, am->lambda / am->mu, am->subsolverZ->solution));
4286285c0a3SHansol  Suh       break;
4296285c0a3SHansol  Suh     }
4306285c0a3SHansol  Suh     am->Hzbool = PETSC_TRUE;
4316285c0a3SHansol  Suh     /* Returns Ax + Bz - c with updated Ax,Bz vectors */
4329566063dSJacob Faibussowitsch     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
4336285c0a3SHansol  Suh     /* Dual variable, y += y + mu*(Ax+Bz-c) */
4349566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(am->y, am->mu, am->residual, am->yold));
4356285c0a3SHansol  Suh 
4366285c0a3SHansol  Suh     /* stopping tolerance update */
4379566063dSJacob Faibussowitsch     PetscCall(TaoADMMToleranceUpdate(tao));
4386285c0a3SHansol  Suh 
4396285c0a3SHansol  Suh     /* Updating Spectral Penalty */
4406285c0a3SHansol  Suh     switch (am->update) {
441d71ae5a4SJacob Faibussowitsch     case TAO_ADMM_UPDATE_BASIC:
442d71ae5a4SJacob Faibussowitsch       am->muold = am->mu;
443d71ae5a4SJacob Faibussowitsch       break;
4446285c0a3SHansol  Suh     case TAO_ADMM_UPDATE_ADAPTIVE:
4456285c0a3SHansol  Suh     case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
4466285c0a3SHansol  Suh       if (tao->niter == 0) {
4479566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->y, am->y0));
4489566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
4491baa6e33SBarry Smith         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
4509566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold));
4519566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Ax, am->Axold));
4529566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Bz, am->Bz0));
4536285c0a3SHansol  Suh         am->muold = am->mu;
4546285c0a3SHansol  Suh       } else if (tao->niter % am->T == 1) {
4556285c0a3SHansol  Suh         /* we have compute Bzold in a previous iteration, and we computed Ax above */
4569566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
4571baa6e33SBarry Smith         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
4589566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(am->yhat, -am->mu, am->residual, am->yold));
4599566063dSJacob Faibussowitsch         PetscCall(AdaptiveADMMPenaltyUpdate(tao));
4609566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Ax, am->Axold));
4619566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->Bz, am->Bz0));
4629566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->yhat, am->yhatold));
4639566063dSJacob Faibussowitsch         PetscCall(VecCopy(am->y, am->y0));
4646285c0a3SHansol  Suh       } else {
4656285c0a3SHansol  Suh         am->muold = am->mu;
4666285c0a3SHansol  Suh       }
4676285c0a3SHansol  Suh       break;
468d71ae5a4SJacob Faibussowitsch     default:
469d71ae5a4SJacob Faibussowitsch       break;
4706285c0a3SHansol  Suh     }
4716285c0a3SHansol  Suh     tao->niter++;
4726285c0a3SHansol  Suh 
4736285c0a3SHansol  Suh     /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/
4746285c0a3SHansol  Suh     switch (am->regswitch) {
4756285c0a3SHansol  Suh     case TAO_ADMM_REGULARIZER_USER:
4766285c0a3SHansol  Suh       if (is_reg_shell) {
4779566063dSJacob Faibussowitsch         PetscCall(ADMML1EpsilonNorm(tao, am->subsolverZ->solution, am->l1epsilon, &reg_func));
4786285c0a3SHansol  Suh       } else {
4793ba16761SJacob Faibussowitsch         PetscCall((*am->ops->regobjgrad)(am->subsolverZ, am->subsolverX->solution, &reg_func, tempL, am->regobjgradP));
4806285c0a3SHansol  Suh       }
4816285c0a3SHansol  Suh       break;
482d71ae5a4SJacob Faibussowitsch     case TAO_ADMM_REGULARIZER_SOFT_THRESH:
483d71ae5a4SJacob Faibussowitsch       PetscCall(ADMML1EpsilonNorm(tao, am->subsolverZ->solution, am->l1epsilon, &reg_func));
484d71ae5a4SJacob Faibussowitsch       break;
4856285c0a3SHansol  Suh     }
4869566063dSJacob Faibussowitsch     PetscCall(VecCopy(am->y, am->yold));
4879566063dSJacob Faibussowitsch     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
4889566063dSJacob Faibussowitsch     PetscCall(VecNorm(am->residual, NORM_2, &am->resnorm));
4899566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, am->last_misfit_val + reg_func, am->dualres, am->resnorm, tao->ksp_its));
4906285c0a3SHansol  Suh 
4919566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, am->last_misfit_val + reg_func, am->dualres, am->resnorm, 1.0));
492dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
4936285c0a3SHansol  Suh   }
4946285c0a3SHansol  Suh   /* Update vectors */
4959566063dSJacob Faibussowitsch   PetscCall(VecCopy(am->subsolverX->solution, tao->solution));
4969566063dSJacob Faibussowitsch   PetscCall(VecCopy(am->subsolverX->gradient, tao->gradient));
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX, "TaoGetADMMParentTao_ADMM", NULL));
4989566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ, "TaoGetADMMParentTao_ADMM", NULL));
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", NULL));
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", NULL));
5019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", NULL));
5029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", NULL));
5033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5046285c0a3SHansol  Suh }
5056285c0a3SHansol  Suh 
506d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_ADMM(Tao tao, PetscOptionItems *PetscOptionsObject)
507d71ae5a4SJacob Faibussowitsch {
5086285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
5096285c0a3SHansol  Suh 
5106285c0a3SHansol  Suh   PetscFunctionBegin;
511d0609cedSBarry 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. ");
5129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient", "regularizer constant", "", am->lambda, &am->lambda, NULL));
5139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty", "Constant for Augmented Lagrangian term.", "", am->mu, &am->mu, NULL));
5149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter", "x relaxation parameter for Z update.", "", am->gamma, &am->gamma, NULL));
5159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor", "ADMM dynamic tolerance update factor.", "", am->tol, &am->tol, NULL));
5169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor", "ADMM spectral penalty update curvature safeguard value.", "", am->orthval, &am->orthval, NULL));
5179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty", "Set ADMM minimum spectral penalty.", "", am->mumin, &am->mumin, NULL));
518d0609cedSBarry Smith   PetscCall(PetscOptionsEnum("-tao_admm_dual_update", "Lagrangian dual update policy", "TaoADMMUpdateType", TaoADMMUpdateTypes, (PetscEnum)am->update, (PetscEnum *)&am->update, NULL));
519d0609cedSBarry Smith   PetscCall(PetscOptionsEnum("-tao_admm_regularizer_type", "ADMM regularizer update rule", "TaoADMMRegularizerType", TaoADMMRegularizerTypes, (PetscEnum)am->regswitch, (PetscEnum *)&am->regswitch, NULL));
520d0609cedSBarry Smith   PetscOptionsHeadEnd();
5219566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(am->subsolverX));
52248a46eb9SPierre Jolivet   if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) PetscCall(TaoSetFromOptions(am->subsolverZ));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5246285c0a3SHansol  Suh }
5256285c0a3SHansol  Suh 
526d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_ADMM(Tao tao, PetscViewer viewer)
527d71ae5a4SJacob Faibussowitsch {
5286285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
5296285c0a3SHansol  Suh 
5306285c0a3SHansol  Suh   PetscFunctionBegin;
5319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
5329566063dSJacob Faibussowitsch   PetscCall(TaoView(am->subsolverX, viewer));
5339566063dSJacob Faibussowitsch   PetscCall(TaoView(am->subsolverZ, viewer));
5349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
5353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5366285c0a3SHansol  Suh }
5376285c0a3SHansol  Suh 
538d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_ADMM(Tao tao)
539d71ae5a4SJacob Faibussowitsch {
5406285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
5416285c0a3SHansol  Suh   PetscInt  n, N, M;
5426285c0a3SHansol  Suh 
5436285c0a3SHansol  Suh   PetscFunctionBegin;
5449566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution, &n));
5459566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &N));
5466285c0a3SHansol  Suh   /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */
5476285c0a3SHansol  Suh   if (!am->JB) {
5486285c0a3SHansol  Suh     am->zJI = PETSC_TRUE;
5499566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao), n, n, PETSC_DETERMINE, PETSC_DETERMINE, NULL, &am->JB));
5509566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JB, MATOP_MULT, (void (*)(void))JacobianIdentityB));
5519566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JB, MATOP_MULT_TRANSPOSE, (void (*)(void))JacobianIdentityB));
5526285c0a3SHansol  Suh     am->JBpre = am->JB;
5536285c0a3SHansol  Suh   }
5546285c0a3SHansol  Suh   if (!am->JA) {
5556285c0a3SHansol  Suh     am->xJI = PETSC_TRUE;
5569566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao), n, n, PETSC_DETERMINE, PETSC_DETERMINE, NULL, &am->JA));
5579566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JA, MATOP_MULT, (void (*)(void))JacobianIdentity));
5589566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(am->JA, MATOP_MULT_TRANSPOSE, (void (*)(void))JacobianIdentity));
5596285c0a3SHansol  Suh     am->JApre = am->JA;
5606285c0a3SHansol  Suh   }
5619566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(am->JA, NULL, &am->Ax));
56248a46eb9SPierre Jolivet   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
5639566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(am->subsolverX, tao->solution));
5646285c0a3SHansol  Suh   if (!am->z) {
5659566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &am->z));
5669566063dSJacob Faibussowitsch     PetscCall(VecSet(am->z, 0.0));
5676285c0a3SHansol  Suh   }
5689566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(am->subsolverZ, am->z));
56948a46eb9SPierre Jolivet   if (!am->workLeft) PetscCall(VecDuplicate(tao->solution, &am->workLeft));
57048a46eb9SPierre Jolivet   if (!am->Axold) PetscCall(VecDuplicate(am->Ax, &am->Axold));
57148a46eb9SPierre Jolivet   if (!am->workJacobianRight) PetscCall(VecDuplicate(am->Ax, &am->workJacobianRight));
57248a46eb9SPierre Jolivet   if (!am->workJacobianRight2) PetscCall(VecDuplicate(am->Ax, &am->workJacobianRight2));
57348a46eb9SPierre Jolivet   if (!am->Bz) PetscCall(VecDuplicate(am->Ax, &am->Bz));
57448a46eb9SPierre Jolivet   if (!am->Bzold) PetscCall(VecDuplicate(am->Ax, &am->Bzold));
57548a46eb9SPierre Jolivet   if (!am->Bz0) PetscCall(VecDuplicate(am->Ax, &am->Bz0));
5766285c0a3SHansol  Suh   if (!am->y) {
5779566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->y));
5789566063dSJacob Faibussowitsch     PetscCall(VecSet(am->y, 0.0));
5796285c0a3SHansol  Suh   }
5806285c0a3SHansol  Suh   if (!am->yold) {
5819566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->yold));
5829566063dSJacob Faibussowitsch     PetscCall(VecSet(am->yold, 0.0));
5836285c0a3SHansol  Suh   }
5846285c0a3SHansol  Suh   if (!am->y0) {
5859566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->y0));
5869566063dSJacob Faibussowitsch     PetscCall(VecSet(am->y0, 0.0));
5876285c0a3SHansol  Suh   }
5886285c0a3SHansol  Suh   if (!am->yhat) {
5899566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->yhat));
5909566063dSJacob Faibussowitsch     PetscCall(VecSet(am->yhat, 0.0));
5916285c0a3SHansol  Suh   }
5926285c0a3SHansol  Suh   if (!am->yhatold) {
5939566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->yhatold));
5949566063dSJacob Faibussowitsch     PetscCall(VecSet(am->yhatold, 0.0));
5956285c0a3SHansol  Suh   }
5966285c0a3SHansol  Suh   if (!am->residual) {
5979566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(am->Ax, &am->residual));
5989566063dSJacob Faibussowitsch     PetscCall(VecSet(am->residual, 0.0));
5996285c0a3SHansol  Suh   }
6006285c0a3SHansol  Suh   if (!am->constraint) {
6016285c0a3SHansol  Suh     am->constraint = NULL;
6026285c0a3SHansol  Suh   } else {
6039566063dSJacob Faibussowitsch     PetscCall(VecGetSize(am->constraint, &M));
6043c859ba3SBarry Smith     PetscCheck(M == N, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Solution vector and constraint vector must be of same size!");
6056285c0a3SHansol  Suh   }
6066285c0a3SHansol  Suh 
6076285c0a3SHansol  Suh   /* Save changed tao tolerance for adaptive tolerance */
608ad540459SPierre Jolivet   if (tao->gatol_changed) am->gatol_admm = tao->gatol;
609ad540459SPierre Jolivet   if (tao->catol_changed) am->catol_admm = tao->catol;
6106285c0a3SHansol  Suh 
6116285c0a3SHansol  Suh   /*Update spectral and dual elements to X subsolver */
6129566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao));
6139566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX, am->JA, am->JApre, am->ops->misfitjac, am->misfitjacobianP));
6149566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ, am->JB, am->JBpre, am->ops->regjac, am->regjacobianP));
6153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6166285c0a3SHansol  Suh }
6176285c0a3SHansol  Suh 
618d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_ADMM(Tao tao)
619d71ae5a4SJacob Faibussowitsch {
6206285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
6216285c0a3SHansol  Suh 
6226285c0a3SHansol  Suh   PetscFunctionBegin;
6239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->z));
6249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Ax));
6259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Axold));
6269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Bz));
6279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Bzold));
6289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->Bz0));
6299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->residual));
6309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->y));
6319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->yold));
6329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->y0));
6339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->yhat));
6349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->yhatold));
6359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->workLeft));
6369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->workJacobianRight));
6379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&am->workJacobianRight2));
6386285c0a3SHansol  Suh 
6399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->JA));
6409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->JB));
64148a46eb9SPierre Jolivet   if (!am->xJI) PetscCall(MatDestroy(&am->JApre));
64248a46eb9SPierre Jolivet   if (!am->zJI) PetscCall(MatDestroy(&am->JBpre));
6436285c0a3SHansol  Suh   if (am->Hx) {
6449566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hx));
6459566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hxpre));
6466285c0a3SHansol  Suh   }
6476285c0a3SHansol  Suh   if (am->Hz) {
6489566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hz));
6499566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hzpre));
6506285c0a3SHansol  Suh   }
6519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->ATA));
6529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&am->BTB));
6539566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&am->subsolverX));
6549566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&am->subsolverZ));
6556285c0a3SHansol  Suh   am->parent = NULL;
6562e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", NULL));
6572e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", NULL));
6582e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", NULL));
6592e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", NULL));
6609566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6626285c0a3SHansol  Suh }
6636285c0a3SHansol  Suh 
6646285c0a3SHansol  Suh /*MC
6656285c0a3SHansol  Suh 
6666285c0a3SHansol  Suh   TAOADMM - Alternating direction method of multipliers method fo solving linear problems with
6676285c0a3SHansol  Suh             constraints. in a min_x f(x) + g(z)  s.t. Ax+Bz=c.
668a5b23f4aSJose E. Roman             This algorithm employs two sub Tao solvers, of which type can be specified
6696285c0a3SHansol  Suh             by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers.
6706285c0a3SHansol  Suh             Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via
6716285c0a3SHansol  Suh             TaoADMMSet{Misfit,Regularizer}HessianChangeStatus.
6726285c0a3SHansol  Suh             Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type.
6736285c0a3SHansol  Suh             There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update,
674a5b23f4aSJose E. Roman             currently there are basic option and adaptive option.
6756285c0a3SHansol  Suh             Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian.
6766285c0a3SHansol  Suh             c can be set with TaoADMMSetConstraintVectorRHS.
6776285c0a3SHansol  Suh             The user can also provide regularizer weight for second subsolver.
6786285c0a3SHansol  Suh 
6796285c0a3SHansol  Suh   References:
680606c0280SSatish Balay . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom
6816285c0a3SHansol  Suh           "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation"
6826285c0a3SHansol  Suh           The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017.
6836285c0a3SHansol  Suh 
6846285c0a3SHansol  Suh   Options Database Keys:
6856285c0a3SHansol  Suh + -tao_admm_regularizer_coefficient        - regularizer constant (default 1.e-6)
6866285c0a3SHansol  Suh . -tao_admm_spectral_penalty               - Constant for Augmented Lagrangian term (default 1.)
6876285c0a3SHansol  Suh . -tao_admm_relaxation_parameter           - relaxation parameter for Z update (default 1.)
6886285c0a3SHansol  Suh . -tao_admm_tolerance_update_factor        - ADMM dynamic tolerance update factor (default 1.e-12)
6896285c0a3SHansol  Suh . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2)
6906285c0a3SHansol  Suh . -tao_admm_minimum_spectral_penalty       - Set ADMM minimum spectral penalty (default 0)
6916285c0a3SHansol  Suh . -tao_admm_dual_update                    - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic")
6926285c0a3SHansol  Suh - -tao_admm_regularizer_type               - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold")
6936285c0a3SHansol  Suh 
6946285c0a3SHansol  Suh   Level: beginner
6956285c0a3SHansol  Suh 
696db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitHessianChangeStatus()`, `TaoADMMSetRegHessianChangeStatus()`, `TaoADMMGetSpectralPenalty()`,
697db781477SPatrick Sanan           `TaoADMMGetMisfitSubsolver()`, `TaoADMMGetRegularizationSubsolver()`, `TaoADMMSetConstraintVectorRHS()`,
698*b4623fecSHansol Suh           `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetRegularizerCoefficient()`, `TaoADMMGetRegularizerCoefficient()`,
699db781477SPatrick Sanan           `TaoADMMSetRegularizerConstraintJacobian()`, `TaoADMMSetMisfitConstraintJacobian()`,
700db781477SPatrick Sanan           `TaoADMMSetMisfitObjectiveAndGradientRoutine()`, `TaoADMMSetMisfitHessianRoutine()`,
701db781477SPatrick Sanan           `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, `TaoADMMSetRegularizerHessianRoutine()`,
702db781477SPatrick Sanan           `TaoGetADMMParentTao()`, `TaoADMMGetDualVector()`, `TaoADMMSetRegularizerType()`,
703db781477SPatrick Sanan           `TaoADMMGetRegularizerType()`, `TaoADMMSetUpdateType()`, `TaoADMMGetUpdateType()`
7046285c0a3SHansol  Suh M*/
7056285c0a3SHansol  Suh 
706d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao)
707d71ae5a4SJacob Faibussowitsch {
7086285c0a3SHansol  Suh   TAO_ADMM *am;
7096285c0a3SHansol  Suh 
7106285c0a3SHansol  Suh   PetscFunctionBegin;
7114dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&am));
7126285c0a3SHansol  Suh 
7136285c0a3SHansol  Suh   tao->ops->destroy        = TaoDestroy_ADMM;
7146285c0a3SHansol  Suh   tao->ops->setup          = TaoSetUp_ADMM;
7156285c0a3SHansol  Suh   tao->ops->setfromoptions = TaoSetFromOptions_ADMM;
7166285c0a3SHansol  Suh   tao->ops->view           = TaoView_ADMM;
7176285c0a3SHansol  Suh   tao->ops->solve          = TaoSolve_ADMM;
7186285c0a3SHansol  Suh 
7196285c0a3SHansol  Suh   tao->data           = (void *)am;
7206285c0a3SHansol  Suh   am->l1epsilon       = 1e-6;
7216285c0a3SHansol  Suh   am->lambda          = 1e-4;
7226285c0a3SHansol  Suh   am->mu              = 1.;
7236285c0a3SHansol  Suh   am->muold           = 0.;
7246285c0a3SHansol  Suh   am->mueps           = PETSC_MACHINE_EPSILON;
7256285c0a3SHansol  Suh   am->mumin           = 0.;
7266285c0a3SHansol  Suh   am->orthval         = 0.2;
7276285c0a3SHansol  Suh   am->T               = 2;
7286285c0a3SHansol  Suh   am->parent          = tao;
7296285c0a3SHansol  Suh   am->update          = TAO_ADMM_UPDATE_BASIC;
7306285c0a3SHansol  Suh   am->regswitch       = TAO_ADMM_REGULARIZER_SOFT_THRESH;
7316285c0a3SHansol  Suh   am->tol             = PETSC_SMALL;
7326285c0a3SHansol  Suh   am->const_norm      = 0;
7336285c0a3SHansol  Suh   am->resnorm         = 0;
7346285c0a3SHansol  Suh   am->dualres         = 0;
73583c8fe1dSLisandro Dalcin   am->ops->regobjgrad = NULL;
73683c8fe1dSLisandro Dalcin   am->ops->reghess    = NULL;
7376285c0a3SHansol  Suh   am->gamma           = 1;
73883c8fe1dSLisandro Dalcin   am->regobjgradP     = NULL;
73983c8fe1dSLisandro Dalcin   am->reghessP        = NULL;
7406285c0a3SHansol  Suh   am->gatol_admm      = 1e-8;
7416285c0a3SHansol  Suh   am->catol_admm      = 0;
7426285c0a3SHansol  Suh   am->Hxchange        = PETSC_TRUE;
7436285c0a3SHansol  Suh   am->Hzchange        = PETSC_TRUE;
7446285c0a3SHansol  Suh   am->Hzbool          = PETSC_TRUE;
7456285c0a3SHansol  Suh   am->Hxbool          = PETSC_TRUE;
7466285c0a3SHansol  Suh 
7479566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &am->subsolverX));
7489566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(am->subsolverX, "misfit_"));
7499566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX, (PetscObject)tao, 1));
7509566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &am->subsolverZ));
7519566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(am->subsolverZ, "reg_"));
7529566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ, (PetscObject)tao, 1));
7536285c0a3SHansol  Suh 
7549566063dSJacob Faibussowitsch   PetscCall(TaoSetType(am->subsolverX, TAONLS));
7559566063dSJacob Faibussowitsch   PetscCall(TaoSetType(am->subsolverZ, TAONLS));
7569566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX, "TaoGetADMMParentTao_ADMM", (PetscObject)tao));
7579566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ, "TaoGetADMMParentTao_ADMM", (PetscObject)tao));
7589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", TaoADMMSetRegularizerType_ADMM));
7599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", TaoADMMGetRegularizerType_ADMM));
7609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", TaoADMMSetUpdateType_ADMM));
7619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", TaoADMMGetUpdateType_ADMM));
7623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7636285c0a3SHansol  Suh }
7646285c0a3SHansol  Suh 
7656285c0a3SHansol  Suh /*@
7666285c0a3SHansol  Suh   TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines  whether Hessian matrix of misfit subsolver changes with respect to input vector.
7676285c0a3SHansol  Suh 
768c3339decSBarry Smith   Collective
7696285c0a3SHansol  Suh 
7706285c0a3SHansol  Suh   Input Parameters:
7717f5c9be9SBarry Smith +  tao - the Tao solver context.
7722fe279fdSBarry Smith -  b - the Hessian matrix change status boolean, `PETSC_FALSE`  when the Hessian matrix does not change, `PETSC_TRUE` otherwise.
7736285c0a3SHansol  Suh 
7746285c0a3SHansol  Suh   Level: advanced
775a5a2f7acSBarry Smith 
776db781477SPatrick Sanan .seealso: `TAOADMM`
7776285c0a3SHansol  Suh @*/
778d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b)
779d71ae5a4SJacob Faibussowitsch {
7806285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
7816285c0a3SHansol  Suh 
7826285c0a3SHansol  Suh   PetscFunctionBegin;
7836285c0a3SHansol  Suh   am->Hxchange = b;
7843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7856285c0a3SHansol  Suh }
7866285c0a3SHansol  Suh 
7876285c0a3SHansol  Suh /*@
7886285c0a3SHansol  Suh   TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector.
7897f5c9be9SBarry Smith 
790c3339decSBarry Smith   Collective
7916285c0a3SHansol  Suh 
7926285c0a3SHansol  Suh   Input Parameters:
7932fe279fdSBarry Smith +  tao - the `Tao` solver context
7942fe279fdSBarry Smith -  b - the Hessian matrix change status boolean, `PETSC_FALSE` when the Hessian matrix does not change, `PETSC_TRUE` otherwise.
7956285c0a3SHansol  Suh 
7966285c0a3SHansol  Suh   Level: advanced
797a5a2f7acSBarry Smith 
798db781477SPatrick Sanan .seealso: `TAOADMM`
7996285c0a3SHansol  Suh @*/
800d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b)
801d71ae5a4SJacob Faibussowitsch {
8026285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8036285c0a3SHansol  Suh 
8046285c0a3SHansol  Suh   PetscFunctionBegin;
8056285c0a3SHansol  Suh   am->Hzchange = b;
8063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8076285c0a3SHansol  Suh }
8086285c0a3SHansol  Suh 
8096285c0a3SHansol  Suh /*@
8106285c0a3SHansol  Suh   TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value
8116285c0a3SHansol  Suh 
812c3339decSBarry Smith   Collective
8136285c0a3SHansol  Suh 
8146285c0a3SHansol  Suh   Input Parameters:
8152fe279fdSBarry Smith +  tao - the `Tao` solver context
8167f5c9be9SBarry Smith -  mu - spectral penalty
8176285c0a3SHansol  Suh 
8186285c0a3SHansol  Suh   Level: advanced
8196285c0a3SHansol  Suh 
820db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TAOADMM`
8216285c0a3SHansol  Suh @*/
822d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu)
823d71ae5a4SJacob Faibussowitsch {
8246285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8256285c0a3SHansol  Suh 
8266285c0a3SHansol  Suh   PetscFunctionBegin;
8276285c0a3SHansol  Suh   am->mu = mu;
8283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8296285c0a3SHansol  Suh }
8306285c0a3SHansol  Suh 
8316285c0a3SHansol  Suh /*@
8326285c0a3SHansol  Suh   TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value
8336285c0a3SHansol  Suh 
834c3339decSBarry Smith   Collective
8356285c0a3SHansol  Suh 
8367f5c9be9SBarry Smith   Input Parameter:
8372fe279fdSBarry Smith .  tao - the `Tao` solver context
8387f5c9be9SBarry Smith 
8397f5c9be9SBarry Smith   Output Parameter:
8407f5c9be9SBarry Smith .  mu - spectral penalty
8416285c0a3SHansol  Suh 
8426285c0a3SHansol  Suh   Level: advanced
8436285c0a3SHansol  Suh 
844db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetSpectralPenalty()`, `TAOADMM`
8456285c0a3SHansol  Suh @*/
846d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu)
847d71ae5a4SJacob Faibussowitsch {
8486285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8496285c0a3SHansol  Suh 
8506285c0a3SHansol  Suh   PetscFunctionBegin;
8516285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
8526285c0a3SHansol  Suh   PetscValidRealPointer(mu, 2);
8536285c0a3SHansol  Suh   *mu = am->mu;
8543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8556285c0a3SHansol  Suh }
8566285c0a3SHansol  Suh 
8576285c0a3SHansol  Suh /*@
8582fe279fdSBarry Smith   TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside `TAOADMM`
8596285c0a3SHansol  Suh 
860c3339decSBarry Smith   Collective
8616285c0a3SHansol  Suh 
8627f5c9be9SBarry Smith   Input Parameter:
8632fe279fdSBarry Smith .  tao - the `Tao` solver context
8647f5c9be9SBarry Smith 
8657f5c9be9SBarry Smith    Output Parameter:
8662fe279fdSBarry Smith .  misfit - the `Tao` subsolver context
8676285c0a3SHansol  Suh 
8686285c0a3SHansol  Suh   Level: advanced
869a5a2f7acSBarry Smith 
8702fe279fdSBarry Smith .seealso: `TAOADMM`, `Tao`
8716285c0a3SHansol  Suh @*/
872d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit)
873d71ae5a4SJacob Faibussowitsch {
8746285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8756285c0a3SHansol  Suh 
8766285c0a3SHansol  Suh   PetscFunctionBegin;
8776285c0a3SHansol  Suh   *misfit = am->subsolverX;
8783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8796285c0a3SHansol  Suh }
8806285c0a3SHansol  Suh 
8816285c0a3SHansol  Suh /*@
8822fe279fdSBarry Smith   TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside `TAOADMM`
8836285c0a3SHansol  Suh 
884c3339decSBarry Smith   Collective
8856285c0a3SHansol  Suh 
8867f5c9be9SBarry Smith   Input Parameter:
8872fe279fdSBarry Smith .  tao - the `Tao` solver context
8887f5c9be9SBarry Smith 
8897f5c9be9SBarry Smith   Output Parameter:
8902fe279fdSBarry Smith .  reg - the `Tao` subsolver context
8916285c0a3SHansol  Suh 
8926285c0a3SHansol  Suh   Level: advanced
893a5a2f7acSBarry Smith 
8942fe279fdSBarry Smith .seealso: `TAOADMM`, `Tao`
8956285c0a3SHansol  Suh @*/
896d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg)
897d71ae5a4SJacob Faibussowitsch {
8986285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
8996285c0a3SHansol  Suh 
9006285c0a3SHansol  Suh   PetscFunctionBegin;
9016285c0a3SHansol  Suh   *reg = am->subsolverZ;
9023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9036285c0a3SHansol  Suh }
9046285c0a3SHansol  Suh 
9056285c0a3SHansol  Suh /*@
9062fe279fdSBarry Smith   TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for `TAOADMM`
9076285c0a3SHansol  Suh 
908c3339decSBarry Smith   Collective
9096285c0a3SHansol  Suh 
9106285c0a3SHansol  Suh   Input Parameters:
9112fe279fdSBarry Smith + tao - the `Tao` solver context
9126285c0a3SHansol  Suh - c - RHS vector
9136285c0a3SHansol  Suh 
9146285c0a3SHansol  Suh   Level: advanced
915a5a2f7acSBarry Smith 
916db781477SPatrick Sanan .seealso: `TAOADMM`
9176285c0a3SHansol  Suh @*/
918d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c)
919d71ae5a4SJacob Faibussowitsch {
9206285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
9216285c0a3SHansol  Suh 
9226285c0a3SHansol  Suh   PetscFunctionBegin;
9236285c0a3SHansol  Suh   am->constraint = c;
9243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9256285c0a3SHansol  Suh }
9266285c0a3SHansol  Suh 
9276285c0a3SHansol  Suh /*@
9287f5c9be9SBarry Smith   TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty
9296285c0a3SHansol  Suh 
930c3339decSBarry Smith   Collective
9316285c0a3SHansol  Suh 
9326285c0a3SHansol  Suh   Input Parameters:
9332fe279fdSBarry Smith +  tao - the `Tao` solver context
9346285c0a3SHansol  Suh -  mu  - minimum spectral penalty value
9356285c0a3SHansol  Suh 
9366285c0a3SHansol  Suh   Level: advanced
9376285c0a3SHansol  Suh 
938db781477SPatrick Sanan .seealso: `TaoADMMGetSpectralPenalty()`, `TAOADMM`
9396285c0a3SHansol  Suh @*/
940d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu)
941d71ae5a4SJacob Faibussowitsch {
9426285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
9436285c0a3SHansol  Suh 
9446285c0a3SHansol  Suh   PetscFunctionBegin;
9456285c0a3SHansol  Suh   am->mumin = mu;
9463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9476285c0a3SHansol  Suh }
9486285c0a3SHansol  Suh 
9496285c0a3SHansol  Suh /*@
9506285c0a3SHansol  Suh   TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case
9516285c0a3SHansol  Suh 
952c3339decSBarry Smith   Collective
9536285c0a3SHansol  Suh 
9546285c0a3SHansol  Suh   Input Parameters:
9552fe279fdSBarry Smith +  tao - the `Tao` solver context
9566285c0a3SHansol  Suh -  lambda - L1-norm regularizer coefficient
9576285c0a3SHansol  Suh 
9586285c0a3SHansol  Suh   Level: advanced
9597f5c9be9SBarry Smith 
960db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
9616285c0a3SHansol  Suh @*/
962d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda)
963d71ae5a4SJacob Faibussowitsch {
9646285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
9656285c0a3SHansol  Suh 
9666285c0a3SHansol  Suh   PetscFunctionBegin;
9676285c0a3SHansol  Suh   am->lambda = lambda;
9683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9696285c0a3SHansol  Suh }
9706285c0a3SHansol  Suh 
971*b4623fecSHansol Suh /*@
972*b4623fecSHansol Suh   TaoADMMGetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case
973*b4623fecSHansol Suh 
974*b4623fecSHansol Suh   Collective
975*b4623fecSHansol Suh 
976*b4623fecSHansol Suh   Input Parameters:
977*b4623fecSHansol Suh +  tao - the `Tao` solver context
978*b4623fecSHansol Suh -  lambda - L1-norm regularizer coefficient
979*b4623fecSHansol Suh 
980*b4623fecSHansol Suh   Level: advanced
981*b4623fecSHansol Suh 
982*b4623fecSHansol Suh .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
983*b4623fecSHansol Suh @*/
984*b4623fecSHansol Suh PetscErrorCode TaoADMMGetRegularizerCoefficient(Tao tao, PetscReal *lambda)
985*b4623fecSHansol Suh {
986*b4623fecSHansol Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
987*b4623fecSHansol Suh 
988*b4623fecSHansol Suh   PetscFunctionBegin;
989*b4623fecSHansol Suh   *lambda = am->lambda;
990*b4623fecSHansol Suh   PetscFunctionReturn(PETSC_SUCCESS);
991*b4623fecSHansol Suh }
992*b4623fecSHansol Suh 
9936285c0a3SHansol  Suh /*@C
9942fe279fdSBarry Smith   TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the `TAOADMM` algorithm. Matrix B constrains the z variable.
9956285c0a3SHansol  Suh 
996c3339decSBarry Smith   Collective
9976285c0a3SHansol  Suh 
998a5a2f7acSBarry Smith   Input Parameters:
999a5a2f7acSBarry Smith + tao - the Tao solver context
10006285c0a3SHansol  Suh . J - user-created regularizer constraint Jacobian matrix
10012fe279fdSBarry Smith . Jpre - user-created regularizer Jacobian constraint matrix for constructing the preconditioner, often this is `J`
10026285c0a3SHansol  Suh . func - function pointer for the regularizer constraint Jacobian update function
10036285c0a3SHansol  Suh - ctx - user context for the regularizer Hessian
10046285c0a3SHansol  Suh 
10056285c0a3SHansol  Suh   Level: advanced
10067f5c9be9SBarry Smith 
1007db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
10086285c0a3SHansol  Suh @*/
1009d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx)
1010d71ae5a4SJacob Faibussowitsch {
10116285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
10126285c0a3SHansol  Suh 
10136285c0a3SHansol  Suh   PetscFunctionBegin;
10146285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10156285c0a3SHansol  Suh   if (J) {
10166285c0a3SHansol  Suh     PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
10176285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, J, 2);
10186285c0a3SHansol  Suh   }
10196285c0a3SHansol  Suh   if (Jpre) {
10206285c0a3SHansol  Suh     PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
10216285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Jpre, 3);
10226285c0a3SHansol  Suh   }
10236285c0a3SHansol  Suh   if (ctx) am->misfitjacobianP = ctx;
10246285c0a3SHansol  Suh   if (func) am->ops->misfitjac = func;
10256285c0a3SHansol  Suh 
10266285c0a3SHansol  Suh   if (J) {
10279566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)J));
10289566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JA));
10296285c0a3SHansol  Suh     am->JA = J;
10306285c0a3SHansol  Suh   }
10316285c0a3SHansol  Suh   if (Jpre) {
10329566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Jpre));
10339566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JApre));
10346285c0a3SHansol  Suh     am->JApre = Jpre;
10356285c0a3SHansol  Suh   }
10363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10376285c0a3SHansol  Suh }
10386285c0a3SHansol  Suh 
10396285c0a3SHansol  Suh /*@C
10402fe279fdSBarry Smith   TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for `TAOADMM` algorithm. Matrix B constraints z variable.
10416285c0a3SHansol  Suh 
1042c3339decSBarry Smith   Collective
10436285c0a3SHansol  Suh 
10446285c0a3SHansol  Suh   Input Parameters:
10452fe279fdSBarry Smith + tao - the `Tao` solver context
10466285c0a3SHansol  Suh . J - user-created regularizer constraint Jacobian matrix
10472fe279fdSBarry Smith . Jpre - user-created regularizer Jacobian constraint matrix for constructing the preconditioner, often this is `J`
10486285c0a3SHansol  Suh . func - function pointer for the regularizer constraint Jacobian update function
10496285c0a3SHansol  Suh - ctx - user context for the regularizer Hessian
10506285c0a3SHansol  Suh 
10516285c0a3SHansol  Suh   Level: advanced
10527f5c9be9SBarry Smith 
1053db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetMisfitConstraintJacobian()`, `TAOADMM`
10546285c0a3SHansol  Suh @*/
1055d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx)
1056d71ae5a4SJacob Faibussowitsch {
10576285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
10586285c0a3SHansol  Suh 
10596285c0a3SHansol  Suh   PetscFunctionBegin;
10606285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10616285c0a3SHansol  Suh   if (J) {
10626285c0a3SHansol  Suh     PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
10636285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, J, 2);
10646285c0a3SHansol  Suh   }
10656285c0a3SHansol  Suh   if (Jpre) {
10666285c0a3SHansol  Suh     PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
10676285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Jpre, 3);
10686285c0a3SHansol  Suh   }
10696285c0a3SHansol  Suh   if (ctx) am->regjacobianP = ctx;
10706285c0a3SHansol  Suh   if (func) am->ops->regjac = func;
10716285c0a3SHansol  Suh 
10726285c0a3SHansol  Suh   if (J) {
10739566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)J));
10749566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JB));
10756285c0a3SHansol  Suh     am->JB = J;
10766285c0a3SHansol  Suh   }
10776285c0a3SHansol  Suh   if (Jpre) {
10789566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Jpre));
10799566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->JBpre));
10806285c0a3SHansol  Suh     am->JBpre = Jpre;
10816285c0a3SHansol  Suh   }
10823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10836285c0a3SHansol  Suh }
10846285c0a3SHansol  Suh 
10856285c0a3SHansol  Suh /*@C
10867f5c9be9SBarry Smith    TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function
10877f5c9be9SBarry Smith 
1088c3339decSBarry Smith    Collective
10896285c0a3SHansol  Suh 
10906285c0a3SHansol  Suh    Input Parameters:
10912fe279fdSBarry Smith +    tao - the `Tao` context
10926285c0a3SHansol  Suh .    func - function pointer for the misfit value and gradient evaluation
10936285c0a3SHansol  Suh -    ctx - user context for the misfit
10946285c0a3SHansol  Suh 
10956285c0a3SHansol  Suh    Level: advanced
10967f5c9be9SBarry Smith 
1097db781477SPatrick Sanan .seealso: `TAOADMM`
10986285c0a3SHansol  Suh @*/
1099d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), void *ctx)
1100d71ae5a4SJacob Faibussowitsch {
11016285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
11026285c0a3SHansol  Suh 
11036285c0a3SHansol  Suh   PetscFunctionBegin;
11046285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11056285c0a3SHansol  Suh   am->misfitobjgradP     = ctx;
11066285c0a3SHansol  Suh   am->ops->misfitobjgrad = func;
11073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11086285c0a3SHansol  Suh }
11096285c0a3SHansol  Suh 
11106285c0a3SHansol  Suh /*@C
11116285c0a3SHansol  Suh    TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back
11126285c0a3SHansol  Suh    function into the algorithm, to be used for subsolverX.
11136285c0a3SHansol  Suh 
1114c3339decSBarry Smith    Collective
11157f5c9be9SBarry Smith 
11166285c0a3SHansol  Suh    Input Parameters:
11172fe279fdSBarry Smith + tao - the `Tao` context
11186285c0a3SHansol  Suh . H - user-created matrix for the Hessian of the misfit term
11196285c0a3SHansol  Suh . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term
11206285c0a3SHansol  Suh . func - function pointer for the misfit Hessian evaluation
11216285c0a3SHansol  Suh - ctx - user context for the misfit Hessian
11226285c0a3SHansol  Suh 
11236285c0a3SHansol  Suh    Level: advanced
1124a5a2f7acSBarry Smith 
1125db781477SPatrick Sanan .seealso: `TAOADMM`
11266285c0a3SHansol  Suh @*/
1127d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx)
1128d71ae5a4SJacob Faibussowitsch {
11296285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
11306285c0a3SHansol  Suh 
11316285c0a3SHansol  Suh   PetscFunctionBegin;
11326285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11336285c0a3SHansol  Suh   if (H) {
11346285c0a3SHansol  Suh     PetscValidHeaderSpecific(H, MAT_CLASSID, 2);
11356285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, H, 2);
11366285c0a3SHansol  Suh   }
11376285c0a3SHansol  Suh   if (Hpre) {
11386285c0a3SHansol  Suh     PetscValidHeaderSpecific(Hpre, MAT_CLASSID, 3);
11396285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Hpre, 3);
11406285c0a3SHansol  Suh   }
1141ad540459SPierre Jolivet   if (ctx) am->misfithessP = ctx;
1142ad540459SPierre Jolivet   if (func) am->ops->misfithess = func;
11436285c0a3SHansol  Suh   if (H) {
11449566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)H));
11459566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hx));
11466285c0a3SHansol  Suh     am->Hx = H;
11476285c0a3SHansol  Suh   }
11486285c0a3SHansol  Suh   if (Hpre) {
11499566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Hpre));
11509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hxpre));
11516285c0a3SHansol  Suh     am->Hxpre = Hpre;
11526285c0a3SHansol  Suh   }
11533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11546285c0a3SHansol  Suh }
11556285c0a3SHansol  Suh 
11566285c0a3SHansol  Suh /*@C
11577f5c9be9SBarry Smith    TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function
11587f5c9be9SBarry Smith 
1159c3339decSBarry Smith    Collective
11606285c0a3SHansol  Suh 
11616285c0a3SHansol  Suh    Input Parameters:
11626285c0a3SHansol  Suh + tao - the Tao context
11636285c0a3SHansol  Suh . func - function pointer for the regularizer value and gradient evaluation
11646285c0a3SHansol  Suh - ctx - user context for the regularizer
11656285c0a3SHansol  Suh 
11666285c0a3SHansol  Suh    Level: advanced
1167a5a2f7acSBarry Smith 
1168db781477SPatrick Sanan .seealso: `TAOADMM`
11696285c0a3SHansol  Suh @*/
1170d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), void *ctx)
1171d71ae5a4SJacob Faibussowitsch {
11726285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
11736285c0a3SHansol  Suh 
11746285c0a3SHansol  Suh   PetscFunctionBegin;
11756285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11766285c0a3SHansol  Suh   am->regobjgradP     = ctx;
11776285c0a3SHansol  Suh   am->ops->regobjgrad = func;
11783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11796285c0a3SHansol  Suh }
11806285c0a3SHansol  Suh 
11816285c0a3SHansol  Suh /*@C
11826285c0a3SHansol  Suh    TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back
11837f5c9be9SBarry Smith    function, to be used for subsolverZ.
11847f5c9be9SBarry Smith 
1185c3339decSBarry Smith    Collective
11866285c0a3SHansol  Suh 
11876285c0a3SHansol  Suh    Input Parameters:
11882fe279fdSBarry Smith + tao - the `Tao` context
11896285c0a3SHansol  Suh . H - user-created matrix for the Hessian of the regularization term
11906285c0a3SHansol  Suh . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term
11916285c0a3SHansol  Suh . func - function pointer for the regularizer Hessian evaluation
11926285c0a3SHansol  Suh - ctx - user context for the regularizer Hessian
11936285c0a3SHansol  Suh 
11946285c0a3SHansol  Suh    Level: advanced
1195a5a2f7acSBarry Smith 
1196db781477SPatrick Sanan .seealso: `TAOADMM`
11976285c0a3SHansol  Suh @*/
1198d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), void *ctx)
1199d71ae5a4SJacob Faibussowitsch {
12006285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
12016285c0a3SHansol  Suh 
12026285c0a3SHansol  Suh   PetscFunctionBegin;
12036285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12046285c0a3SHansol  Suh   if (H) {
12056285c0a3SHansol  Suh     PetscValidHeaderSpecific(H, MAT_CLASSID, 2);
12066285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, H, 2);
12076285c0a3SHansol  Suh   }
12086285c0a3SHansol  Suh   if (Hpre) {
12096285c0a3SHansol  Suh     PetscValidHeaderSpecific(Hpre, MAT_CLASSID, 3);
12106285c0a3SHansol  Suh     PetscCheckSameComm(tao, 1, Hpre, 3);
12116285c0a3SHansol  Suh   }
1212ad540459SPierre Jolivet   if (ctx) am->reghessP = ctx;
1213ad540459SPierre Jolivet   if (func) am->ops->reghess = func;
12146285c0a3SHansol  Suh   if (H) {
12159566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)H));
12169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hz));
12176285c0a3SHansol  Suh     am->Hz = H;
12186285c0a3SHansol  Suh   }
12196285c0a3SHansol  Suh   if (Hpre) {
12209566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Hpre));
12219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&am->Hzpre));
12226285c0a3SHansol  Suh     am->Hzpre = Hpre;
12236285c0a3SHansol  Suh   }
12243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12256285c0a3SHansol  Suh }
12266285c0a3SHansol  Suh 
12276285c0a3SHansol  Suh /*@
12282fe279fdSBarry Smith    TaoGetADMMParentTao - Gets pointer to parent `TAOADMM`, used by inner subsolver.
12296285c0a3SHansol  Suh 
1230c3339decSBarry Smith    Collective
12317f5c9be9SBarry Smith 
12327f5c9be9SBarry Smith    Input Parameter:
12332fe279fdSBarry Smith . tao - the `Tao` context
12347f5c9be9SBarry Smith 
12357f5c9be9SBarry Smith    Output Parameter:
12362fe279fdSBarry Smith . admm_tao - the parent `Tao` context
12376285c0a3SHansol  Suh 
12386285c0a3SHansol  Suh    Level: advanced
1239a5a2f7acSBarry Smith 
1240db781477SPatrick Sanan .seealso: `TAOADMM`
12416285c0a3SHansol  Suh @*/
1242d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao)
1243d71ae5a4SJacob Faibussowitsch {
12446285c0a3SHansol  Suh   PetscFunctionBegin;
12456285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12469566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)tao, "TaoGetADMMParentTao_ADMM", (PetscObject *)admm_tao));
12473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12486285c0a3SHansol  Suh }
12496285c0a3SHansol  Suh 
12506285c0a3SHansol  Suh /*@
12512fe279fdSBarry Smith   TaoADMMGetDualVector - Returns the dual vector associated with the current `TAOADMM` state
12526285c0a3SHansol  Suh 
12536285c0a3SHansol  Suh   Not Collective
12546285c0a3SHansol  Suh 
12556285c0a3SHansol  Suh   Input Parameter:
12562fe279fdSBarry Smith . tao - the `Tao` context
12577f5c9be9SBarry Smith 
12587f5c9be9SBarry Smith   Output Parameter:
12597f5c9be9SBarry Smith . Y - the current solution
12606285c0a3SHansol  Suh 
12616285c0a3SHansol  Suh   Level: intermediate
1262a5a2f7acSBarry Smith 
1263db781477SPatrick Sanan .seealso: `TAOADMM`
12646285c0a3SHansol  Suh @*/
1265d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y)
1266d71ae5a4SJacob Faibussowitsch {
12676285c0a3SHansol  Suh   TAO_ADMM *am = (TAO_ADMM *)tao->data;
12686285c0a3SHansol  Suh 
12696285c0a3SHansol  Suh   PetscFunctionBegin;
12706285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12716285c0a3SHansol  Suh   *Y = am->y;
12723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12736285c0a3SHansol  Suh }
12746285c0a3SHansol  Suh 
12756285c0a3SHansol  Suh /*@
12762fe279fdSBarry Smith   TaoADMMSetRegularizerType - Set regularizer type for `TAOADMM` routine
12776285c0a3SHansol  Suh 
12786285c0a3SHansol  Suh   Not Collective
12796285c0a3SHansol  Suh 
1280d8d19677SJose E. Roman   Input Parameters:
12812fe279fdSBarry Smith + tao  - the `Tao` context
12827f5c9be9SBarry Smith - type - regularizer type
12837f5c9be9SBarry Smith 
12843c7db156SBarry Smith   Options Database Key:
1285147403d9SBarry Smith .  -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer
12866285c0a3SHansol  Suh 
12876285c0a3SHansol  Suh   Level: intermediate
12886285c0a3SHansol  Suh 
1289db781477SPatrick Sanan .seealso: `TaoADMMGetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
12906285c0a3SHansol  Suh @*/
1291d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type)
1292d71ae5a4SJacob Faibussowitsch {
12936285c0a3SHansol  Suh   PetscFunctionBegin;
12946285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12956285c0a3SHansol  Suh   PetscValidLogicalCollectiveEnum(tao, type, 2);
1296cac4c232SBarry Smith   PetscTryMethod(tao, "TaoADMMSetRegularizerType_C", (Tao, TaoADMMRegularizerType), (tao, type));
12973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12986285c0a3SHansol  Suh }
12996285c0a3SHansol  Suh 
13006285c0a3SHansol  Suh /*@
13012fe279fdSBarry Smith    TaoADMMGetRegularizerType - Gets the type of regularizer routine for `TAOADMM`
13026285c0a3SHansol  Suh 
13036285c0a3SHansol  Suh    Not Collective
13046285c0a3SHansol  Suh 
13056285c0a3SHansol  Suh    Input Parameter:
13062fe279fdSBarry Smith .  tao - the `Tao` context
13076285c0a3SHansol  Suh 
13086285c0a3SHansol  Suh    Output Parameter:
13096285c0a3SHansol  Suh .  type - the type of regularizer
13106285c0a3SHansol  Suh 
13116285c0a3SHansol  Suh    Level: intermediate
13126285c0a3SHansol  Suh 
1313db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
13146285c0a3SHansol  Suh @*/
1315d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type)
1316d71ae5a4SJacob Faibussowitsch {
13176285c0a3SHansol  Suh   PetscFunctionBegin;
13186285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1319cac4c232SBarry Smith   PetscUseMethod(tao, "TaoADMMGetRegularizerType_C", (Tao, TaoADMMRegularizerType *), (tao, type));
13203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13216285c0a3SHansol  Suh }
13226285c0a3SHansol  Suh 
13236285c0a3SHansol  Suh /*@
13242fe279fdSBarry Smith   TaoADMMSetUpdateType - Set update routine for `TAOADMM` routine
13256285c0a3SHansol  Suh 
13266285c0a3SHansol  Suh   Not Collective
13276285c0a3SHansol  Suh 
1328d8d19677SJose E. Roman   Input Parameters:
13292fe279fdSBarry Smith + tao  - the `Tao` context
13306285c0a3SHansol  Suh - type - spectral parameter update type
13316285c0a3SHansol  Suh 
13326285c0a3SHansol  Suh   Level: intermediate
13336285c0a3SHansol  Suh 
1334db781477SPatrick Sanan .seealso: `TaoADMMGetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
13356285c0a3SHansol  Suh @*/
1336d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type)
1337d71ae5a4SJacob Faibussowitsch {
13386285c0a3SHansol  Suh   PetscFunctionBegin;
13396285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
13406285c0a3SHansol  Suh   PetscValidLogicalCollectiveEnum(tao, type, 2);
1341cac4c232SBarry Smith   PetscTryMethod(tao, "TaoADMMSetUpdateType_C", (Tao, TaoADMMUpdateType), (tao, type));
13423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13436285c0a3SHansol  Suh }
13446285c0a3SHansol  Suh 
13456285c0a3SHansol  Suh /*@
13462fe279fdSBarry Smith    TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for `TAOADMM`
13476285c0a3SHansol  Suh 
13486285c0a3SHansol  Suh    Not Collective
13496285c0a3SHansol  Suh 
13506285c0a3SHansol  Suh    Input Parameter:
13512fe279fdSBarry Smith .  tao - the `Tao` context
13526285c0a3SHansol  Suh 
13536285c0a3SHansol  Suh    Output Parameter:
13546285c0a3SHansol  Suh .  type - the type of spectral penalty update routine
13556285c0a3SHansol  Suh 
13566285c0a3SHansol  Suh    Level: intermediate
13576285c0a3SHansol  Suh 
1358db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
13596285c0a3SHansol  Suh @*/
1360d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type)
1361d71ae5a4SJacob Faibussowitsch {
13626285c0a3SHansol  Suh   PetscFunctionBegin;
13636285c0a3SHansol  Suh   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1364cac4c232SBarry Smith   PetscUseMethod(tao, "TaoADMMGetUpdateType_C", (Tao, TaoADMMUpdateType *), (tao, type));
13653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13666285c0a3SHansol  Suh }
1367