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
TaoADMMToleranceUpdate(Tao tao)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));
56606f75f6SBarry Smith PetscCall(TaoSetConstraintTolerances(tao, temp, PETSC_CURRENT));
57606f75f6SBarry Smith PetscCall(TaoSetTolerances(tao, am->gatol_admm * ATynorm, PETSC_CURRENT, PETSC_CURRENT));
583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
596285c0a3SHansol Suh }
606285c0a3SHansol Suh
6146091a0eSPierre Jolivet /* Penalty Update for Adaptive ADMM. */
AdaptiveADMMPenaltyUpdate(Tao tao)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
TaoADMMSetRegularizerType_ADMM(Tao tao,TaoADMMRegularizerType type)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
TaoADMMGetRegularizerType_ADMM(Tao tao,TaoADMMRegularizerType * type)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
TaoADMMSetUpdateType_ADMM(Tao tao,TaoADMMUpdateType type)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
TaoADMMGetUpdateType_ADMM(Tao tao,TaoADMMUpdateType * type)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*/
ADMMUpdateConstraintResidualVector(Tao tao,Vec x,Vec z,Vec Ax,Vec Bz,Vec residual)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. */
SubObjGradUpdate(Tao tao,Vec x,PetscReal * f,Vec g,void * ptr)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. */
RegObjGradUpdate(Tao tao,Vec z,PetscReal * f,Vec g,void * ptr)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 */
ADMML1EpsilonNorm(Tao tao,Vec x,PetscReal eps,PetscReal * norm)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
ADMMInternalHessianUpdate(Mat H,Mat Constraint,PetscBool Identity,void * ptr)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) {
250f971d498SPierre Jolivet case TAO_ADMM_UPDATE_BASIC:
251d71ae5a4SJacob Faibussowitsch break;
252f971d498SPierre Jolivet case TAO_ADMM_UPDATE_ADAPTIVE:
253f971d498SPierre Jolivet 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 */
SubHessianUpdate(Tao tao,Vec x,Mat H,Mat Hpre,void * ptr)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 */
28515229ffcSPierre Jolivet /* Update Lagrangian only once 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 */
RegHessianUpdate(Tao tao,Vec z,Mat H,Mat Hpre,void * ptr)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 if (am->Hzchange) {
3006285c0a3SHansol Suh /* Case where Hessian gets updated with respect to x vector input. */
3019566063dSJacob Faibussowitsch PetscCall((*am->ops->reghess)(am->subsolverZ, z, H, Hpre, am->reghessP));
3029566063dSJacob Faibussowitsch PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian, am->BTB, am->zJI, am));
3036285c0a3SHansol Suh } else if (am->Hzbool) {
3046285c0a3SHansol Suh /* Hessian doesn't get updated. H(x) = c */
30515229ffcSPierre Jolivet /* Update Lagrangian only once per TAO call */
3069566063dSJacob Faibussowitsch PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian, am->BTB, am->zJI, am));
3076285c0a3SHansol Suh am->Hzbool = PETSC_FALSE;
3086285c0a3SHansol Suh }
3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3106285c0a3SHansol Suh }
3116285c0a3SHansol Suh
3126285c0a3SHansol Suh /* Shell Matrix routine for A matrix.
3136285c0a3SHansol Suh * This gets used when user puts NULL for
3146285c0a3SHansol Suh * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...)
3156285c0a3SHansol Suh * Essentially sets A=I*/
JacobianIdentity(Mat mat,Vec in,Vec out)316d71ae5a4SJacob Faibussowitsch static PetscErrorCode JacobianIdentity(Mat mat, Vec in, Vec out)
317d71ae5a4SJacob Faibussowitsch {
3186285c0a3SHansol Suh PetscFunctionBegin;
3199566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out));
3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3216285c0a3SHansol Suh }
3226285c0a3SHansol Suh
3236285c0a3SHansol Suh /* Shell Matrix routine for B matrix.
3246285c0a3SHansol Suh * This gets used when user puts NULL for
3256285c0a3SHansol Suh * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...)
3266285c0a3SHansol Suh * Sets B=-I */
JacobianIdentityB(Mat mat,Vec in,Vec out)327d71ae5a4SJacob Faibussowitsch static PetscErrorCode JacobianIdentityB(Mat mat, Vec in, Vec out)
328d71ae5a4SJacob Faibussowitsch {
3296285c0a3SHansol Suh PetscFunctionBegin;
3309566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out));
3319566063dSJacob Faibussowitsch PetscCall(VecScale(out, -1.));
3323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3336285c0a3SHansol Suh }
3346285c0a3SHansol Suh
3356285c0a3SHansol Suh /* Solve f(x) + g(z) s.t. Ax + Bz = c */
TaoSolve_ADMM(Tao tao)336d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_ADMM(Tao tao)
337d71ae5a4SJacob Faibussowitsch {
3386285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
3396285c0a3SHansol Suh PetscInt N;
3406285c0a3SHansol Suh PetscReal reg_func;
3416285c0a3SHansol Suh PetscBool is_reg_shell;
3426285c0a3SHansol Suh Vec tempL;
3436285c0a3SHansol Suh
3446285c0a3SHansol Suh PetscFunctionBegin;
3456285c0a3SHansol Suh if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) {
3463c859ba3SBarry Smith PetscCheck(am->subsolverX->ops->computejacobianequality, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoADMMSetMisfitConstraintJacobian() first");
3473c859ba3SBarry Smith PetscCheck(am->subsolverZ->ops->computejacobianequality, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoADMMSetRegularizerConstraintJacobian() first");
34848a46eb9SPierre Jolivet if (am->constraint != NULL) PetscCall(VecNorm(am->constraint, NORM_2, &am->const_norm));
3496285c0a3SHansol Suh }
3506285c0a3SHansol Suh tempL = am->workLeft;
3519566063dSJacob Faibussowitsch PetscCall(VecGetSize(tempL, &N));
3526285c0a3SHansol Suh
35348a46eb9SPierre Jolivet if (am->Hx && am->ops->misfithess) PetscCall(TaoSetHessian(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao));
3546285c0a3SHansol Suh
3556285c0a3SHansol Suh if (!am->zJI) {
3566285c0a3SHansol Suh /* Currently, B is assumed to be a linear system, i.e., not getting updated*/
357fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(am->JB, am->JB, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &am->BTB));
3586285c0a3SHansol Suh }
3596285c0a3SHansol Suh if (!am->xJI) {
3606285c0a3SHansol Suh /* Currently, A is assumed to be a linear system, i.e., not getting updated*/
361fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &am->ATA));
3626285c0a3SHansol Suh }
3636285c0a3SHansol Suh
3646285c0a3SHansol Suh is_reg_shell = PETSC_FALSE;
3656285c0a3SHansol Suh
3669566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell));
3676285c0a3SHansol Suh
3686285c0a3SHansol Suh if (!is_reg_shell) {
3696285c0a3SHansol Suh switch (am->regswitch) {
370f971d498SPierre Jolivet case TAO_ADMM_REGULARIZER_USER:
371d71ae5a4SJacob Faibussowitsch break;
372f971d498SPierre Jolivet case TAO_ADMM_REGULARIZER_SOFT_THRESH:
3736285c0a3SHansol Suh /* Soft Threshold. */
3746285c0a3SHansol Suh break;
3756285c0a3SHansol Suh }
3761baa6e33SBarry Smith if (am->ops->regobjgrad) PetscCall(TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao));
37748a46eb9SPierre Jolivet if (am->Hz && am->ops->reghess) PetscCall(TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao));
3786285c0a3SHansol Suh }
3796285c0a3SHansol Suh
3806285c0a3SHansol Suh switch (am->update) {
3816285c0a3SHansol Suh case TAO_ADMM_UPDATE_BASIC:
3826285c0a3SHansol Suh if (am->subsolverX->hessian) {
3836285c0a3SHansol Suh /* In basic case, Hessian does not get updated w.r.t. to spectral penalty
3846285c0a3SHansol Suh * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/
3856285c0a3SHansol Suh if (!am->xJI) {
3869566063dSJacob Faibussowitsch PetscCall(MatAXPY(am->subsolverX->hessian, am->mu, am->ATA, DIFFERENT_NONZERO_PATTERN));
3876285c0a3SHansol Suh } else {
3889566063dSJacob Faibussowitsch PetscCall(MatShift(am->subsolverX->hessian, am->mu));
3896285c0a3SHansol Suh }
3906285c0a3SHansol Suh }
3916285c0a3SHansol Suh if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) {
3926285c0a3SHansol Suh if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) {
3939566063dSJacob Faibussowitsch PetscCall(MatAXPY(am->subsolverZ->hessian, am->mu, am->BTB, DIFFERENT_NONZERO_PATTERN));
3946285c0a3SHansol Suh } else {
3959566063dSJacob Faibussowitsch PetscCall(MatShift(am->subsolverZ->hessian, am->mu));
3966285c0a3SHansol Suh }
3976285c0a3SHansol Suh }
3986285c0a3SHansol Suh break;
3996285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE:
400d71ae5a4SJacob Faibussowitsch case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
401d71ae5a4SJacob Faibussowitsch break;
4026285c0a3SHansol Suh }
4036285c0a3SHansol Suh
4049566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation, &cited));
4056285c0a3SHansol Suh tao->reason = TAO_CONTINUE_ITERATING;
4066285c0a3SHansol Suh
4076285c0a3SHansol Suh while (tao->reason == TAO_CONTINUE_ITERATING) {
408dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
4099566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Bz, am->Bzold));
4106285c0a3SHansol Suh
4116285c0a3SHansol Suh /* x update */
4129566063dSJacob Faibussowitsch PetscCall(TaoSolve(am->subsolverX));
4139566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre));
4149566063dSJacob Faibussowitsch PetscCall(MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution, am->Ax));
4156285c0a3SHansol Suh
4166285c0a3SHansol Suh am->Hxbool = PETSC_TRUE;
4176285c0a3SHansol Suh
4186285c0a3SHansol Suh /* z update */
4196285c0a3SHansol Suh switch (am->regswitch) {
420d71ae5a4SJacob Faibussowitsch case TAO_ADMM_REGULARIZER_USER:
421d71ae5a4SJacob Faibussowitsch PetscCall(TaoSolve(am->subsolverZ));
422d71ae5a4SJacob Faibussowitsch break;
4236285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_SOFT_THRESH:
4246285c0a3SHansol Suh /* L1 assumes A,B jacobians are identity nxn matrix */
4259566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->workJacobianRight, 1 / am->mu, am->y, am->Ax));
4269566063dSJacob Faibussowitsch PetscCall(TaoSoftThreshold(am->workJacobianRight, -am->lambda / am->mu, am->lambda / am->mu, am->subsolverZ->solution));
4276285c0a3SHansol Suh break;
4286285c0a3SHansol Suh }
4296285c0a3SHansol Suh am->Hzbool = PETSC_TRUE;
4306285c0a3SHansol Suh /* Returns Ax + Bz - c with updated Ax,Bz vectors */
4319566063dSJacob Faibussowitsch PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
4326285c0a3SHansol Suh /* Dual variable, y += y + mu*(Ax+Bz-c) */
4339566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->y, am->mu, am->residual, am->yold));
4346285c0a3SHansol Suh
4356285c0a3SHansol Suh /* stopping tolerance update */
4369566063dSJacob Faibussowitsch PetscCall(TaoADMMToleranceUpdate(tao));
4376285c0a3SHansol Suh
4386285c0a3SHansol Suh /* Updating Spectral Penalty */
4396285c0a3SHansol Suh switch (am->update) {
440d71ae5a4SJacob Faibussowitsch case TAO_ADMM_UPDATE_BASIC:
441d71ae5a4SJacob Faibussowitsch am->muold = am->mu;
442d71ae5a4SJacob Faibussowitsch break;
4436285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE:
4446285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
4456285c0a3SHansol Suh if (tao->niter == 0) {
4469566063dSJacob Faibussowitsch PetscCall(VecCopy(am->y, am->y0));
4479566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
4481baa6e33SBarry Smith if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
4499566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold));
4509566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Ax, am->Axold));
4519566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Bz, am->Bz0));
4526285c0a3SHansol Suh am->muold = am->mu;
4536285c0a3SHansol Suh } else if (tao->niter % am->T == 1) {
4546285c0a3SHansol Suh /* we have compute Bzold in a previous iteration, and we computed Ax above */
4559566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
4561baa6e33SBarry Smith if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
4579566063dSJacob Faibussowitsch PetscCall(VecWAXPY(am->yhat, -am->mu, am->residual, am->yold));
4589566063dSJacob Faibussowitsch PetscCall(AdaptiveADMMPenaltyUpdate(tao));
4599566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Ax, am->Axold));
4609566063dSJacob Faibussowitsch PetscCall(VecCopy(am->Bz, am->Bz0));
4619566063dSJacob Faibussowitsch PetscCall(VecCopy(am->yhat, am->yhatold));
4629566063dSJacob Faibussowitsch PetscCall(VecCopy(am->y, am->y0));
4636285c0a3SHansol Suh } else {
4646285c0a3SHansol Suh am->muold = am->mu;
4656285c0a3SHansol Suh }
4666285c0a3SHansol Suh break;
467d71ae5a4SJacob Faibussowitsch default:
468d71ae5a4SJacob Faibussowitsch break;
4696285c0a3SHansol Suh }
4706285c0a3SHansol Suh tao->niter++;
4716285c0a3SHansol Suh
4726285c0a3SHansol Suh /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/
4736285c0a3SHansol Suh switch (am->regswitch) {
4746285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_USER:
4756285c0a3SHansol Suh if (is_reg_shell) {
4769566063dSJacob Faibussowitsch PetscCall(ADMML1EpsilonNorm(tao, am->subsolverZ->solution, am->l1epsilon, ®_func));
4776285c0a3SHansol Suh } else {
4783ba16761SJacob Faibussowitsch PetscCall((*am->ops->regobjgrad)(am->subsolverZ, am->subsolverX->solution, ®_func, tempL, am->regobjgradP));
4796285c0a3SHansol Suh }
4806285c0a3SHansol Suh break;
481d71ae5a4SJacob Faibussowitsch case TAO_ADMM_REGULARIZER_SOFT_THRESH:
482d71ae5a4SJacob Faibussowitsch PetscCall(ADMML1EpsilonNorm(tao, am->subsolverZ->solution, am->l1epsilon, ®_func));
483d71ae5a4SJacob Faibussowitsch break;
4846285c0a3SHansol Suh }
4859566063dSJacob Faibussowitsch PetscCall(VecCopy(am->y, am->yold));
4869566063dSJacob Faibussowitsch PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
4879566063dSJacob Faibussowitsch PetscCall(VecNorm(am->residual, NORM_2, &am->resnorm));
4889566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, am->last_misfit_val + reg_func, am->dualres, am->resnorm, tao->ksp_its));
4896285c0a3SHansol Suh
4909566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, am->last_misfit_val + reg_func, am->dualres, am->resnorm, 1.0));
491dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
4926285c0a3SHansol Suh }
4936285c0a3SHansol Suh /* Update vectors */
4949566063dSJacob Faibussowitsch PetscCall(VecCopy(am->subsolverX->solution, tao->solution));
4959566063dSJacob Faibussowitsch PetscCall(VecCopy(am->subsolverX->gradient, tao->gradient));
4969566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverX, "TaoGetADMMParentTao_ADMM", NULL));
4979566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ, "TaoGetADMMParentTao_ADMM", NULL));
4989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", NULL));
4999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", NULL));
5009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", NULL));
5019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", NULL));
5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5036285c0a3SHansol Suh }
5046285c0a3SHansol Suh
TaoSetFromOptions_ADMM(Tao tao,PetscOptionItems PetscOptionsObject)505ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_ADMM(Tao tao, PetscOptionItems PetscOptionsObject)
506d71ae5a4SJacob Faibussowitsch {
5076285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
5086285c0a3SHansol Suh
5096285c0a3SHansol Suh PetscFunctionBegin;
510d0609cedSBarry 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. ");
5119566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient", "regularizer constant", "", am->lambda, &am->lambda, NULL));
5129566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty", "Constant for Augmented Lagrangian term.", "", am->mu, &am->mu, NULL));
5139566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter", "x relaxation parameter for Z update.", "", am->gamma, &am->gamma, NULL));
5149566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor", "ADMM dynamic tolerance update factor.", "", am->tol, &am->tol, NULL));
5159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor", "ADMM spectral penalty update curvature safeguard value.", "", am->orthval, &am->orthval, NULL));
5169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty", "Set ADMM minimum spectral penalty.", "", am->mumin, &am->mumin, NULL));
517d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-tao_admm_dual_update", "Lagrangian dual update policy", "TaoADMMUpdateType", TaoADMMUpdateTypes, (PetscEnum)am->update, (PetscEnum *)&am->update, NULL));
518d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-tao_admm_regularizer_type", "ADMM regularizer update rule", "TaoADMMRegularizerType", TaoADMMRegularizerTypes, (PetscEnum)am->regswitch, (PetscEnum *)&am->regswitch, NULL));
519d0609cedSBarry Smith PetscOptionsHeadEnd();
5209566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(am->subsolverX));
52148a46eb9SPierre Jolivet if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) PetscCall(TaoSetFromOptions(am->subsolverZ));
5223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5236285c0a3SHansol Suh }
5246285c0a3SHansol Suh
TaoView_ADMM(Tao tao,PetscViewer viewer)525d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_ADMM(Tao tao, PetscViewer viewer)
526d71ae5a4SJacob Faibussowitsch {
5276285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
5286285c0a3SHansol Suh
5296285c0a3SHansol Suh PetscFunctionBegin;
5309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
5319566063dSJacob Faibussowitsch PetscCall(TaoView(am->subsolverX, viewer));
5329566063dSJacob Faibussowitsch PetscCall(TaoView(am->subsolverZ, viewer));
5339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5356285c0a3SHansol Suh }
5366285c0a3SHansol Suh
TaoSetUp_ADMM(Tao tao)537d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_ADMM(Tao tao)
538d71ae5a4SJacob Faibussowitsch {
5396285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
5406285c0a3SHansol Suh PetscInt n, N, M;
5416285c0a3SHansol Suh
5426285c0a3SHansol Suh PetscFunctionBegin;
5439566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n));
5449566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N));
5456285c0a3SHansol Suh /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */
5466285c0a3SHansol Suh if (!am->JB) {
5476285c0a3SHansol Suh am->zJI = PETSC_TRUE;
5489566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao), n, n, PETSC_DETERMINE, PETSC_DETERMINE, NULL, &am->JB));
54957d50842SBarry Smith PetscCall(MatShellSetOperation(am->JB, MATOP_MULT, (PetscErrorCodeFn *)JacobianIdentityB));
55057d50842SBarry Smith PetscCall(MatShellSetOperation(am->JB, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)JacobianIdentityB));
5516285c0a3SHansol Suh am->JBpre = am->JB;
5526285c0a3SHansol Suh }
5536285c0a3SHansol Suh if (!am->JA) {
5546285c0a3SHansol Suh am->xJI = PETSC_TRUE;
5559566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao), n, n, PETSC_DETERMINE, PETSC_DETERMINE, NULL, &am->JA));
55657d50842SBarry Smith PetscCall(MatShellSetOperation(am->JA, MATOP_MULT, (PetscErrorCodeFn *)JacobianIdentity));
55757d50842SBarry Smith PetscCall(MatShellSetOperation(am->JA, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)JacobianIdentity));
5586285c0a3SHansol Suh am->JApre = am->JA;
5596285c0a3SHansol Suh }
5609566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(am->JA, NULL, &am->Ax));
56148a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
5629566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(am->subsolverX, tao->solution));
5636285c0a3SHansol Suh if (!am->z) {
5649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &am->z));
5659566063dSJacob Faibussowitsch PetscCall(VecSet(am->z, 0.0));
5666285c0a3SHansol Suh }
5679566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(am->subsolverZ, am->z));
56848a46eb9SPierre Jolivet if (!am->workLeft) PetscCall(VecDuplicate(tao->solution, &am->workLeft));
56948a46eb9SPierre Jolivet if (!am->Axold) PetscCall(VecDuplicate(am->Ax, &am->Axold));
57048a46eb9SPierre Jolivet if (!am->workJacobianRight) PetscCall(VecDuplicate(am->Ax, &am->workJacobianRight));
57148a46eb9SPierre Jolivet if (!am->workJacobianRight2) PetscCall(VecDuplicate(am->Ax, &am->workJacobianRight2));
57248a46eb9SPierre Jolivet if (!am->Bz) PetscCall(VecDuplicate(am->Ax, &am->Bz));
57348a46eb9SPierre Jolivet if (!am->Bzold) PetscCall(VecDuplicate(am->Ax, &am->Bzold));
57448a46eb9SPierre Jolivet if (!am->Bz0) PetscCall(VecDuplicate(am->Ax, &am->Bz0));
5756285c0a3SHansol Suh if (!am->y) {
5769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax, &am->y));
5779566063dSJacob Faibussowitsch PetscCall(VecSet(am->y, 0.0));
5786285c0a3SHansol Suh }
5796285c0a3SHansol Suh if (!am->yold) {
5809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax, &am->yold));
5819566063dSJacob Faibussowitsch PetscCall(VecSet(am->yold, 0.0));
5826285c0a3SHansol Suh }
5836285c0a3SHansol Suh if (!am->y0) {
5849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax, &am->y0));
5859566063dSJacob Faibussowitsch PetscCall(VecSet(am->y0, 0.0));
5866285c0a3SHansol Suh }
5876285c0a3SHansol Suh if (!am->yhat) {
5889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax, &am->yhat));
5899566063dSJacob Faibussowitsch PetscCall(VecSet(am->yhat, 0.0));
5906285c0a3SHansol Suh }
5916285c0a3SHansol Suh if (!am->yhatold) {
5929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax, &am->yhatold));
5939566063dSJacob Faibussowitsch PetscCall(VecSet(am->yhatold, 0.0));
5946285c0a3SHansol Suh }
5956285c0a3SHansol Suh if (!am->residual) {
5969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(am->Ax, &am->residual));
5979566063dSJacob Faibussowitsch PetscCall(VecSet(am->residual, 0.0));
5986285c0a3SHansol Suh }
5996285c0a3SHansol Suh if (!am->constraint) {
6006285c0a3SHansol Suh am->constraint = NULL;
6016285c0a3SHansol Suh } else {
6029566063dSJacob Faibussowitsch PetscCall(VecGetSize(am->constraint, &M));
6033c859ba3SBarry Smith PetscCheck(M == N, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Solution vector and constraint vector must be of same size!");
6046285c0a3SHansol Suh }
6056285c0a3SHansol Suh
6066285c0a3SHansol Suh /* Save changed tao tolerance for adaptive tolerance */
607606f75f6SBarry Smith if (tao->gatol != tao->default_gatol) am->gatol_admm = tao->gatol;
608606f75f6SBarry Smith if (tao->catol != tao->default_catol) am->catol_admm = tao->catol;
6096285c0a3SHansol Suh
6106285c0a3SHansol Suh /*Update spectral and dual elements to X subsolver */
6119566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao));
6129566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX, am->JA, am->JApre, am->ops->misfitjac, am->misfitjacobianP));
6139566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ, am->JB, am->JBpre, am->ops->regjac, am->regjacobianP));
6143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6156285c0a3SHansol Suh }
6166285c0a3SHansol Suh
TaoDestroy_ADMM(Tao tao)617d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_ADMM(Tao tao)
618d71ae5a4SJacob Faibussowitsch {
6196285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
6206285c0a3SHansol Suh
6216285c0a3SHansol Suh PetscFunctionBegin;
6229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->z));
6239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Ax));
6249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Axold));
6259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Bz));
6269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Bzold));
6279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->Bz0));
6289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->residual));
6299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->y));
6309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->yold));
6319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->y0));
6329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->yhat));
6339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->yhatold));
6349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->workLeft));
6359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->workJacobianRight));
6369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&am->workJacobianRight2));
6376285c0a3SHansol Suh
6389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JA));
6399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JB));
64048a46eb9SPierre Jolivet if (!am->xJI) PetscCall(MatDestroy(&am->JApre));
64148a46eb9SPierre Jolivet if (!am->zJI) PetscCall(MatDestroy(&am->JBpre));
6426285c0a3SHansol Suh if (am->Hx) {
6439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hx));
6449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hxpre));
6456285c0a3SHansol Suh }
6466285c0a3SHansol Suh if (am->Hz) {
6479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hz));
6489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hzpre));
6496285c0a3SHansol Suh }
6509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->ATA));
6519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->BTB));
6529566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&am->subsolverX));
6539566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&am->subsolverZ));
6546285c0a3SHansol Suh am->parent = NULL;
6552e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", NULL));
6562e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", NULL));
6572e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", NULL));
6582e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", NULL));
6599566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
6603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6616285c0a3SHansol Suh }
6626285c0a3SHansol Suh
6636285c0a3SHansol Suh /*MC
6645804573cSPierre Jolivet TAOADMM - Alternating direction method of multipliers method for solving linear problems with
6651d27aa22SBarry Smith constraints. in a $ \min_x f(x) + g(z)$ s.t. $Ax+Bz=c$.
666a5b23f4aSJose E. Roman This algorithm employs two sub Tao solvers, of which type can be specified
6676285c0a3SHansol Suh by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers.
6686285c0a3SHansol Suh Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via
6691d27aa22SBarry Smith `TaoADMMSet{Misfit,Regularizer}HessianChangeStatus()`.
6701d27aa22SBarry Smith Second subsolver does support `TAOSHELL`. It should be noted that L1-norm is used for objective value for `TAOSHELL` type.
6716285c0a3SHansol Suh There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update,
672a5b23f4aSJose E. Roman currently there are basic option and adaptive option.
6731d27aa22SBarry Smith Constraint is set at Ax+Bz=c, and A and B can be set with `TaoADMMSet{Misfit,Regularizer}ConstraintJacobian()`.
6741d27aa22SBarry Smith c can be set with `TaoADMMSetConstraintVectorRHS()`.
6751d27aa22SBarry Smith The user can also provide regularizer weight for second subsolver. {cite}`xu2017adaptive`
6766285c0a3SHansol Suh
6776285c0a3SHansol Suh Options Database Keys:
6786285c0a3SHansol Suh + -tao_admm_regularizer_coefficient - regularizer constant (default 1.e-6)
6796285c0a3SHansol Suh . -tao_admm_spectral_penalty - Constant for Augmented Lagrangian term (default 1.)
6806285c0a3SHansol Suh . -tao_admm_relaxation_parameter - relaxation parameter for Z update (default 1.)
6816285c0a3SHansol Suh . -tao_admm_tolerance_update_factor - ADMM dynamic tolerance update factor (default 1.e-12)
6826285c0a3SHansol Suh . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2)
6836285c0a3SHansol Suh . -tao_admm_minimum_spectral_penalty - Set ADMM minimum spectral penalty (default 0)
6846285c0a3SHansol Suh . -tao_admm_dual_update - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic")
6856285c0a3SHansol Suh - -tao_admm_regularizer_type - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold")
6866285c0a3SHansol Suh
6876285c0a3SHansol Suh Level: beginner
6886285c0a3SHansol Suh
689db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitHessianChangeStatus()`, `TaoADMMSetRegHessianChangeStatus()`, `TaoADMMGetSpectralPenalty()`,
690db781477SPatrick Sanan `TaoADMMGetMisfitSubsolver()`, `TaoADMMGetRegularizationSubsolver()`, `TaoADMMSetConstraintVectorRHS()`,
691b4623fecSHansol Suh `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetRegularizerCoefficient()`, `TaoADMMGetRegularizerCoefficient()`,
692db781477SPatrick Sanan `TaoADMMSetRegularizerConstraintJacobian()`, `TaoADMMSetMisfitConstraintJacobian()`,
693db781477SPatrick Sanan `TaoADMMSetMisfitObjectiveAndGradientRoutine()`, `TaoADMMSetMisfitHessianRoutine()`,
694db781477SPatrick Sanan `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, `TaoADMMSetRegularizerHessianRoutine()`,
695db781477SPatrick Sanan `TaoGetADMMParentTao()`, `TaoADMMGetDualVector()`, `TaoADMMSetRegularizerType()`,
696db781477SPatrick Sanan `TaoADMMGetRegularizerType()`, `TaoADMMSetUpdateType()`, `TaoADMMGetUpdateType()`
6976285c0a3SHansol Suh M*/
6986285c0a3SHansol Suh
TaoCreate_ADMM(Tao tao)699d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao)
700d71ae5a4SJacob Faibussowitsch {
7016285c0a3SHansol Suh TAO_ADMM *am;
7026285c0a3SHansol Suh
7036285c0a3SHansol Suh PetscFunctionBegin;
7044dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&am));
7056285c0a3SHansol Suh
7066285c0a3SHansol Suh tao->ops->destroy = TaoDestroy_ADMM;
7076285c0a3SHansol Suh tao->ops->setup = TaoSetUp_ADMM;
7086285c0a3SHansol Suh tao->ops->setfromoptions = TaoSetFromOptions_ADMM;
7096285c0a3SHansol Suh tao->ops->view = TaoView_ADMM;
7106285c0a3SHansol Suh tao->ops->solve = TaoSolve_ADMM;
7116285c0a3SHansol Suh
712606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
713606f75f6SBarry Smith
7146285c0a3SHansol Suh tao->data = (void *)am;
7156285c0a3SHansol Suh am->l1epsilon = 1e-6;
7166285c0a3SHansol Suh am->lambda = 1e-4;
7176285c0a3SHansol Suh am->mu = 1.;
7186285c0a3SHansol Suh am->muold = 0.;
7196285c0a3SHansol Suh am->mueps = PETSC_MACHINE_EPSILON;
7206285c0a3SHansol Suh am->mumin = 0.;
7216285c0a3SHansol Suh am->orthval = 0.2;
7226285c0a3SHansol Suh am->T = 2;
7236285c0a3SHansol Suh am->parent = tao;
7246285c0a3SHansol Suh am->update = TAO_ADMM_UPDATE_BASIC;
7256285c0a3SHansol Suh am->regswitch = TAO_ADMM_REGULARIZER_SOFT_THRESH;
7266285c0a3SHansol Suh am->tol = PETSC_SMALL;
7276285c0a3SHansol Suh am->const_norm = 0;
7286285c0a3SHansol Suh am->resnorm = 0;
7296285c0a3SHansol Suh am->dualres = 0;
73083c8fe1dSLisandro Dalcin am->ops->regobjgrad = NULL;
73183c8fe1dSLisandro Dalcin am->ops->reghess = NULL;
7326285c0a3SHansol Suh am->gamma = 1;
73383c8fe1dSLisandro Dalcin am->regobjgradP = NULL;
73483c8fe1dSLisandro Dalcin am->reghessP = NULL;
7356285c0a3SHansol Suh am->gatol_admm = 1e-8;
7366285c0a3SHansol Suh am->catol_admm = 0;
7376285c0a3SHansol Suh am->Hxchange = PETSC_TRUE;
7386285c0a3SHansol Suh am->Hzchange = PETSC_TRUE;
7396285c0a3SHansol Suh am->Hzbool = PETSC_TRUE;
7406285c0a3SHansol Suh am->Hxbool = PETSC_TRUE;
7416285c0a3SHansol Suh
7429566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &am->subsolverX));
7439566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(am->subsolverX, "misfit_"));
7449566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX, (PetscObject)tao, 1));
7459566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &am->subsolverZ));
7469566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(am->subsolverZ, "reg_"));
7479566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ, (PetscObject)tao, 1));
7486285c0a3SHansol Suh
7499566063dSJacob Faibussowitsch PetscCall(TaoSetType(am->subsolverX, TAONLS));
7509566063dSJacob Faibussowitsch PetscCall(TaoSetType(am->subsolverZ, TAONLS));
7519566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverX, "TaoGetADMMParentTao_ADMM", (PetscObject)tao));
7529566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ, "TaoGetADMMParentTao_ADMM", (PetscObject)tao));
7539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetRegularizerType_C", TaoADMMSetRegularizerType_ADMM));
7549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetRegularizerType_C", TaoADMMGetRegularizerType_ADMM));
7559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMSetUpdateType_C", TaoADMMSetUpdateType_ADMM));
7569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoADMMGetUpdateType_C", TaoADMMGetUpdateType_ADMM));
7573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7586285c0a3SHansol Suh }
7596285c0a3SHansol Suh
7606285c0a3SHansol Suh /*@
7616285c0a3SHansol Suh TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines whether Hessian matrix of misfit subsolver changes with respect to input vector.
7626285c0a3SHansol Suh
763c3339decSBarry Smith Collective
7646285c0a3SHansol Suh
7656285c0a3SHansol Suh Input Parameters:
7667f5c9be9SBarry Smith + tao - the Tao solver context.
7672fe279fdSBarry Smith - b - the Hessian matrix change status boolean, `PETSC_FALSE` when the Hessian matrix does not change, `PETSC_TRUE` otherwise.
7686285c0a3SHansol Suh
7696285c0a3SHansol Suh Level: advanced
770a5a2f7acSBarry Smith
771db781477SPatrick Sanan .seealso: `TAOADMM`
7726285c0a3SHansol Suh @*/
TaoADMMSetMisfitHessianChangeStatus(Tao tao,PetscBool b)773d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b)
774d71ae5a4SJacob Faibussowitsch {
7756285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
7766285c0a3SHansol Suh
7776285c0a3SHansol Suh PetscFunctionBegin;
7786285c0a3SHansol Suh am->Hxchange = b;
7793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7806285c0a3SHansol Suh }
7816285c0a3SHansol Suh
7826285c0a3SHansol Suh /*@
7836285c0a3SHansol Suh TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector.
7847f5c9be9SBarry Smith
785c3339decSBarry Smith Collective
7866285c0a3SHansol Suh
7876285c0a3SHansol Suh Input Parameters:
7882fe279fdSBarry Smith + tao - the `Tao` solver context
7892fe279fdSBarry Smith - b - the Hessian matrix change status boolean, `PETSC_FALSE` when the Hessian matrix does not change, `PETSC_TRUE` otherwise.
7906285c0a3SHansol Suh
7916285c0a3SHansol Suh Level: advanced
792a5a2f7acSBarry Smith
793db781477SPatrick Sanan .seealso: `TAOADMM`
7946285c0a3SHansol Suh @*/
TaoADMMSetRegHessianChangeStatus(Tao tao,PetscBool b)795d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b)
796d71ae5a4SJacob Faibussowitsch {
7976285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
7986285c0a3SHansol Suh
7996285c0a3SHansol Suh PetscFunctionBegin;
8006285c0a3SHansol Suh am->Hzchange = b;
8013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8026285c0a3SHansol Suh }
8036285c0a3SHansol Suh
8046285c0a3SHansol Suh /*@
8056285c0a3SHansol Suh TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value
8066285c0a3SHansol Suh
807c3339decSBarry Smith Collective
8086285c0a3SHansol Suh
8096285c0a3SHansol Suh Input Parameters:
8102fe279fdSBarry Smith + tao - the `Tao` solver context
8117f5c9be9SBarry Smith - mu - spectral penalty
8126285c0a3SHansol Suh
8136285c0a3SHansol Suh Level: advanced
8146285c0a3SHansol Suh
815db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TAOADMM`
8166285c0a3SHansol Suh @*/
TaoADMMSetSpectralPenalty(Tao tao,PetscReal mu)817d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu)
818d71ae5a4SJacob Faibussowitsch {
8196285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
8206285c0a3SHansol Suh
8216285c0a3SHansol Suh PetscFunctionBegin;
8226285c0a3SHansol Suh am->mu = mu;
8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8246285c0a3SHansol Suh }
8256285c0a3SHansol Suh
8266285c0a3SHansol Suh /*@
8276285c0a3SHansol Suh TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value
8286285c0a3SHansol Suh
829c3339decSBarry Smith Collective
8306285c0a3SHansol Suh
8317f5c9be9SBarry Smith Input Parameter:
8322fe279fdSBarry Smith . tao - the `Tao` solver context
8337f5c9be9SBarry Smith
8347f5c9be9SBarry Smith Output Parameter:
8357f5c9be9SBarry Smith . mu - spectral penalty
8366285c0a3SHansol Suh
8376285c0a3SHansol Suh Level: advanced
8386285c0a3SHansol Suh
839db781477SPatrick Sanan .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetSpectralPenalty()`, `TAOADMM`
8406285c0a3SHansol Suh @*/
TaoADMMGetSpectralPenalty(Tao tao,PetscReal * mu)841d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu)
842d71ae5a4SJacob Faibussowitsch {
8436285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
8446285c0a3SHansol Suh
8456285c0a3SHansol Suh PetscFunctionBegin;
8466285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
8474f572ea9SToby Isaac PetscAssertPointer(mu, 2);
8486285c0a3SHansol Suh *mu = am->mu;
8493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8506285c0a3SHansol Suh }
8516285c0a3SHansol Suh
8526285c0a3SHansol Suh /*@
8532fe279fdSBarry Smith TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside `TAOADMM`
8546285c0a3SHansol Suh
855c3339decSBarry Smith Collective
8566285c0a3SHansol Suh
8577f5c9be9SBarry Smith Input Parameter:
8582fe279fdSBarry Smith . tao - the `Tao` solver context
8597f5c9be9SBarry Smith
8607f5c9be9SBarry Smith Output Parameter:
8612fe279fdSBarry Smith . misfit - the `Tao` subsolver context
8626285c0a3SHansol Suh
8636285c0a3SHansol Suh Level: advanced
864a5a2f7acSBarry Smith
8652fe279fdSBarry Smith .seealso: `TAOADMM`, `Tao`
8666285c0a3SHansol Suh @*/
TaoADMMGetMisfitSubsolver(Tao tao,Tao * misfit)867d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit)
868d71ae5a4SJacob Faibussowitsch {
8696285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
8706285c0a3SHansol Suh
8716285c0a3SHansol Suh PetscFunctionBegin;
8726285c0a3SHansol Suh *misfit = am->subsolverX;
8733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8746285c0a3SHansol Suh }
8756285c0a3SHansol Suh
8766285c0a3SHansol Suh /*@
8772fe279fdSBarry Smith TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside `TAOADMM`
8786285c0a3SHansol Suh
879c3339decSBarry Smith Collective
8806285c0a3SHansol Suh
8817f5c9be9SBarry Smith Input Parameter:
8822fe279fdSBarry Smith . tao - the `Tao` solver context
8837f5c9be9SBarry Smith
8847f5c9be9SBarry Smith Output Parameter:
8852fe279fdSBarry Smith . reg - the `Tao` subsolver context
8866285c0a3SHansol Suh
8876285c0a3SHansol Suh Level: advanced
888a5a2f7acSBarry Smith
8892fe279fdSBarry Smith .seealso: `TAOADMM`, `Tao`
8906285c0a3SHansol Suh @*/
TaoADMMGetRegularizationSubsolver(Tao tao,Tao * reg)891d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg)
892d71ae5a4SJacob Faibussowitsch {
8936285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
8946285c0a3SHansol Suh
8956285c0a3SHansol Suh PetscFunctionBegin;
8966285c0a3SHansol Suh *reg = am->subsolverZ;
8973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8986285c0a3SHansol Suh }
8996285c0a3SHansol Suh
9006285c0a3SHansol Suh /*@
9012fe279fdSBarry Smith TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for `TAOADMM`
9026285c0a3SHansol Suh
903c3339decSBarry Smith Collective
9046285c0a3SHansol Suh
9056285c0a3SHansol Suh Input Parameters:
9062fe279fdSBarry Smith + tao - the `Tao` solver context
9076285c0a3SHansol Suh - c - RHS vector
9086285c0a3SHansol Suh
9096285c0a3SHansol Suh Level: advanced
910a5a2f7acSBarry Smith
911db781477SPatrick Sanan .seealso: `TAOADMM`
9126285c0a3SHansol Suh @*/
TaoADMMSetConstraintVectorRHS(Tao tao,Vec c)913d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c)
914d71ae5a4SJacob Faibussowitsch {
9156285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
9166285c0a3SHansol Suh
9176285c0a3SHansol Suh PetscFunctionBegin;
9186285c0a3SHansol Suh am->constraint = c;
9193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9206285c0a3SHansol Suh }
9216285c0a3SHansol Suh
9226285c0a3SHansol Suh /*@
9237f5c9be9SBarry Smith TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty
9246285c0a3SHansol Suh
925c3339decSBarry Smith Collective
9266285c0a3SHansol Suh
9276285c0a3SHansol Suh Input Parameters:
9282fe279fdSBarry Smith + tao - the `Tao` solver context
9296285c0a3SHansol Suh - mu - minimum spectral penalty value
9306285c0a3SHansol Suh
9316285c0a3SHansol Suh Level: advanced
9326285c0a3SHansol Suh
933db781477SPatrick Sanan .seealso: `TaoADMMGetSpectralPenalty()`, `TAOADMM`
9346285c0a3SHansol Suh @*/
TaoADMMSetMinimumSpectralPenalty(Tao tao,PetscReal mu)935d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu)
936d71ae5a4SJacob Faibussowitsch {
9376285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
9386285c0a3SHansol Suh
9396285c0a3SHansol Suh PetscFunctionBegin;
9406285c0a3SHansol Suh am->mumin = mu;
9413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9426285c0a3SHansol Suh }
9436285c0a3SHansol Suh
9446285c0a3SHansol Suh /*@
9456285c0a3SHansol Suh TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case
9466285c0a3SHansol Suh
947c3339decSBarry Smith Collective
9486285c0a3SHansol Suh
9496285c0a3SHansol Suh Input Parameters:
9502fe279fdSBarry Smith + tao - the `Tao` solver context
9516285c0a3SHansol Suh - lambda - L1-norm regularizer coefficient
9526285c0a3SHansol Suh
9536285c0a3SHansol Suh Level: advanced
9547f5c9be9SBarry Smith
955db781477SPatrick Sanan .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
9566285c0a3SHansol Suh @*/
TaoADMMSetRegularizerCoefficient(Tao tao,PetscReal lambda)957d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda)
958d71ae5a4SJacob Faibussowitsch {
9596285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
9606285c0a3SHansol Suh
9616285c0a3SHansol Suh PetscFunctionBegin;
9626285c0a3SHansol Suh am->lambda = lambda;
9633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9646285c0a3SHansol Suh }
9656285c0a3SHansol Suh
966b4623fecSHansol Suh /*@
967e6c0fc48SBarry Smith TaoADMMGetRegularizerCoefficient - Get the regularization coefficient lambda for L1 norm regularization case
968b4623fecSHansol Suh
969b4623fecSHansol Suh Collective
970b4623fecSHansol Suh
971e6c0fc48SBarry Smith Input Parameter:
972e6c0fc48SBarry Smith . tao - the `Tao` solver context
973e6c0fc48SBarry Smith
974e6c0fc48SBarry Smith Output Parameter:
975e6c0fc48SBarry Smith . lambda - L1-norm regularizer coefficient
976b4623fecSHansol Suh
977b4623fecSHansol Suh Level: advanced
978b4623fecSHansol Suh
979b4623fecSHansol Suh .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
980b4623fecSHansol Suh @*/
TaoADMMGetRegularizerCoefficient(Tao tao,PetscReal * lambda)981b4623fecSHansol Suh PetscErrorCode TaoADMMGetRegularizerCoefficient(Tao tao, PetscReal *lambda)
982b4623fecSHansol Suh {
983b4623fecSHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
984b4623fecSHansol Suh
985b4623fecSHansol Suh PetscFunctionBegin;
986b4623fecSHansol Suh *lambda = am->lambda;
987b4623fecSHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
988b4623fecSHansol Suh }
989b4623fecSHansol Suh
9906285c0a3SHansol Suh /*@C
9912fe279fdSBarry Smith TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the `TAOADMM` algorithm. Matrix B constrains the z variable.
9926285c0a3SHansol Suh
993c3339decSBarry Smith Collective
9946285c0a3SHansol Suh
995a5a2f7acSBarry Smith Input Parameters:
996a5a2f7acSBarry Smith + tao - the Tao solver context
9976285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix
9982fe279fdSBarry Smith . Jpre - user-created regularizer Jacobian constraint matrix for constructing the preconditioner, often this is `J`
9996285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function
10006285c0a3SHansol Suh - ctx - user context for the regularizer Hessian
10016285c0a3SHansol Suh
10026285c0a3SHansol Suh Level: advanced
10037f5c9be9SBarry Smith
1004db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
10056285c0a3SHansol Suh @*/
TaoADMMSetMisfitConstraintJacobian(Tao tao,Mat J,Mat Jpre,PetscErrorCode (* func)(Tao,Vec,Mat,Mat,void *),PetscCtx ctx)1006*2a8381b2SBarry Smith PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), PetscCtx ctx)
1007d71ae5a4SJacob Faibussowitsch {
10086285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
10096285c0a3SHansol Suh
10106285c0a3SHansol Suh PetscFunctionBegin;
10116285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10126285c0a3SHansol Suh if (J) {
10136285c0a3SHansol Suh PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
10146285c0a3SHansol Suh PetscCheckSameComm(tao, 1, J, 2);
10156285c0a3SHansol Suh }
10166285c0a3SHansol Suh if (Jpre) {
10176285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
10186285c0a3SHansol Suh PetscCheckSameComm(tao, 1, Jpre, 3);
10196285c0a3SHansol Suh }
10206285c0a3SHansol Suh if (ctx) am->misfitjacobianP = ctx;
10216285c0a3SHansol Suh if (func) am->ops->misfitjac = func;
10226285c0a3SHansol Suh
10236285c0a3SHansol Suh if (J) {
10249566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
10259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JA));
10266285c0a3SHansol Suh am->JA = J;
10276285c0a3SHansol Suh }
10286285c0a3SHansol Suh if (Jpre) {
10299566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre));
10309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JApre));
10316285c0a3SHansol Suh am->JApre = Jpre;
10326285c0a3SHansol Suh }
10333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10346285c0a3SHansol Suh }
10356285c0a3SHansol Suh
10366285c0a3SHansol Suh /*@C
10372fe279fdSBarry Smith TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for `TAOADMM` algorithm. Matrix B constraints z variable.
10386285c0a3SHansol Suh
1039c3339decSBarry Smith Collective
10406285c0a3SHansol Suh
10416285c0a3SHansol Suh Input Parameters:
10422fe279fdSBarry Smith + tao - the `Tao` solver context
10436285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix
10442fe279fdSBarry Smith . Jpre - user-created regularizer Jacobian constraint matrix for constructing the preconditioner, often this is `J`
10456285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function
10466285c0a3SHansol Suh - ctx - user context for the regularizer Hessian
10476285c0a3SHansol Suh
10486285c0a3SHansol Suh Level: advanced
10497f5c9be9SBarry Smith
1050db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetMisfitConstraintJacobian()`, `TAOADMM`
10516285c0a3SHansol Suh @*/
TaoADMMSetRegularizerConstraintJacobian(Tao tao,Mat J,Mat Jpre,PetscErrorCode (* func)(Tao,Vec,Mat,Mat,void *),PetscCtx ctx)1052*2a8381b2SBarry Smith PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), PetscCtx ctx)
1053d71ae5a4SJacob Faibussowitsch {
10546285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
10556285c0a3SHansol Suh
10566285c0a3SHansol Suh PetscFunctionBegin;
10576285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10586285c0a3SHansol Suh if (J) {
10596285c0a3SHansol Suh PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
10606285c0a3SHansol Suh PetscCheckSameComm(tao, 1, J, 2);
10616285c0a3SHansol Suh }
10626285c0a3SHansol Suh if (Jpre) {
10636285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
10646285c0a3SHansol Suh PetscCheckSameComm(tao, 1, Jpre, 3);
10656285c0a3SHansol Suh }
10666285c0a3SHansol Suh if (ctx) am->regjacobianP = ctx;
10676285c0a3SHansol Suh if (func) am->ops->regjac = func;
10686285c0a3SHansol Suh
10696285c0a3SHansol Suh if (J) {
10709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
10719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JB));
10726285c0a3SHansol Suh am->JB = J;
10736285c0a3SHansol Suh }
10746285c0a3SHansol Suh if (Jpre) {
10759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre));
10769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->JBpre));
10776285c0a3SHansol Suh am->JBpre = Jpre;
10786285c0a3SHansol Suh }
10793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10806285c0a3SHansol Suh }
10816285c0a3SHansol Suh
10826285c0a3SHansol Suh /*@C
10837f5c9be9SBarry Smith TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function
10847f5c9be9SBarry Smith
1085c3339decSBarry Smith Collective
10866285c0a3SHansol Suh
10876285c0a3SHansol Suh Input Parameters:
10882fe279fdSBarry Smith + tao - the `Tao` context
10896285c0a3SHansol Suh . func - function pointer for the misfit value and gradient evaluation
10906285c0a3SHansol Suh - ctx - user context for the misfit
10916285c0a3SHansol Suh
10926285c0a3SHansol Suh Level: advanced
10937f5c9be9SBarry Smith
1094db781477SPatrick Sanan .seealso: `TAOADMM`
10956285c0a3SHansol Suh @*/
TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (* func)(Tao,Vec,PetscReal *,Vec,void *),PetscCtx ctx)1096*2a8381b2SBarry Smith PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), PetscCtx ctx)
1097d71ae5a4SJacob Faibussowitsch {
10986285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
10996285c0a3SHansol Suh
11006285c0a3SHansol Suh PetscFunctionBegin;
11016285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11026285c0a3SHansol Suh am->misfitobjgradP = ctx;
11036285c0a3SHansol Suh am->ops->misfitobjgrad = func;
11043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11056285c0a3SHansol Suh }
11066285c0a3SHansol Suh
11076285c0a3SHansol Suh /*@C
11086285c0a3SHansol Suh TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back
11096285c0a3SHansol Suh function into the algorithm, to be used for subsolverX.
11106285c0a3SHansol Suh
1111c3339decSBarry Smith Collective
11127f5c9be9SBarry Smith
11136285c0a3SHansol Suh Input Parameters:
11142fe279fdSBarry Smith + tao - the `Tao` context
11156285c0a3SHansol Suh . H - user-created matrix for the Hessian of the misfit term
11166285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term
11176285c0a3SHansol Suh . func - function pointer for the misfit Hessian evaluation
11186285c0a3SHansol Suh - ctx - user context for the misfit Hessian
11196285c0a3SHansol Suh
11206285c0a3SHansol Suh Level: advanced
1121a5a2f7acSBarry Smith
1122db781477SPatrick Sanan .seealso: `TAOADMM`
11236285c0a3SHansol Suh @*/
TaoADMMSetMisfitHessianRoutine(Tao tao,Mat H,Mat Hpre,PetscErrorCode (* func)(Tao,Vec,Mat,Mat,void *),PetscCtx ctx)1124*2a8381b2SBarry Smith PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), PetscCtx ctx)
1125d71ae5a4SJacob Faibussowitsch {
11266285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
11276285c0a3SHansol Suh
11286285c0a3SHansol Suh PetscFunctionBegin;
11296285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11306285c0a3SHansol Suh if (H) {
11316285c0a3SHansol Suh PetscValidHeaderSpecific(H, MAT_CLASSID, 2);
11326285c0a3SHansol Suh PetscCheckSameComm(tao, 1, H, 2);
11336285c0a3SHansol Suh }
11346285c0a3SHansol Suh if (Hpre) {
11356285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre, MAT_CLASSID, 3);
11366285c0a3SHansol Suh PetscCheckSameComm(tao, 1, Hpre, 3);
11376285c0a3SHansol Suh }
1138ad540459SPierre Jolivet if (ctx) am->misfithessP = ctx;
1139ad540459SPierre Jolivet if (func) am->ops->misfithess = func;
11406285c0a3SHansol Suh if (H) {
11419566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H));
11429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hx));
11436285c0a3SHansol Suh am->Hx = H;
11446285c0a3SHansol Suh }
11456285c0a3SHansol Suh if (Hpre) {
11469566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Hpre));
11479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hxpre));
11486285c0a3SHansol Suh am->Hxpre = Hpre;
11496285c0a3SHansol Suh }
11503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11516285c0a3SHansol Suh }
11526285c0a3SHansol Suh
11536285c0a3SHansol Suh /*@C
11547f5c9be9SBarry Smith TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function
11557f5c9be9SBarry Smith
1156c3339decSBarry Smith Collective
11576285c0a3SHansol Suh
11586285c0a3SHansol Suh Input Parameters:
11596285c0a3SHansol Suh + tao - the Tao context
11606285c0a3SHansol Suh . func - function pointer for the regularizer value and gradient evaluation
11616285c0a3SHansol Suh - ctx - user context for the regularizer
11626285c0a3SHansol Suh
11636285c0a3SHansol Suh Level: advanced
1164a5a2f7acSBarry Smith
1165db781477SPatrick Sanan .seealso: `TAOADMM`
11666285c0a3SHansol Suh @*/
TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao,PetscErrorCode (* func)(Tao,Vec,PetscReal *,Vec,void *),PetscCtx ctx)1167*2a8381b2SBarry Smith PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), PetscCtx ctx)
1168d71ae5a4SJacob Faibussowitsch {
11696285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
11706285c0a3SHansol Suh
11716285c0a3SHansol Suh PetscFunctionBegin;
11726285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11736285c0a3SHansol Suh am->regobjgradP = ctx;
11746285c0a3SHansol Suh am->ops->regobjgrad = func;
11753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11766285c0a3SHansol Suh }
11776285c0a3SHansol Suh
11786285c0a3SHansol Suh /*@C
11796285c0a3SHansol Suh TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back
11807f5c9be9SBarry Smith function, to be used for subsolverZ.
11817f5c9be9SBarry Smith
1182c3339decSBarry Smith Collective
11836285c0a3SHansol Suh
11846285c0a3SHansol Suh Input Parameters:
11852fe279fdSBarry Smith + tao - the `Tao` context
11866285c0a3SHansol Suh . H - user-created matrix for the Hessian of the regularization term
11876285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term
11886285c0a3SHansol Suh . func - function pointer for the regularizer Hessian evaluation
11896285c0a3SHansol Suh - ctx - user context for the regularizer Hessian
11906285c0a3SHansol Suh
11916285c0a3SHansol Suh Level: advanced
1192a5a2f7acSBarry Smith
1193db781477SPatrick Sanan .seealso: `TAOADMM`
11946285c0a3SHansol Suh @*/
TaoADMMSetRegularizerHessianRoutine(Tao tao,Mat H,Mat Hpre,PetscErrorCode (* func)(Tao,Vec,Mat,Mat,void *),PetscCtx ctx)1195*2a8381b2SBarry Smith PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void *), PetscCtx ctx)
1196d71ae5a4SJacob Faibussowitsch {
11976285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
11986285c0a3SHansol Suh
11996285c0a3SHansol Suh PetscFunctionBegin;
12006285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12016285c0a3SHansol Suh if (H) {
12026285c0a3SHansol Suh PetscValidHeaderSpecific(H, MAT_CLASSID, 2);
12036285c0a3SHansol Suh PetscCheckSameComm(tao, 1, H, 2);
12046285c0a3SHansol Suh }
12056285c0a3SHansol Suh if (Hpre) {
12066285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre, MAT_CLASSID, 3);
12076285c0a3SHansol Suh PetscCheckSameComm(tao, 1, Hpre, 3);
12086285c0a3SHansol Suh }
1209ad540459SPierre Jolivet if (ctx) am->reghessP = ctx;
1210ad540459SPierre Jolivet if (func) am->ops->reghess = func;
12116285c0a3SHansol Suh if (H) {
12129566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H));
12139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hz));
12146285c0a3SHansol Suh am->Hz = H;
12156285c0a3SHansol Suh }
12166285c0a3SHansol Suh if (Hpre) {
12179566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Hpre));
12189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&am->Hzpre));
12196285c0a3SHansol Suh am->Hzpre = Hpre;
12206285c0a3SHansol Suh }
12213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12226285c0a3SHansol Suh }
12236285c0a3SHansol Suh
12246285c0a3SHansol Suh /*@
12252fe279fdSBarry Smith TaoGetADMMParentTao - Gets pointer to parent `TAOADMM`, used by inner subsolver.
12266285c0a3SHansol Suh
1227c3339decSBarry Smith Collective
12287f5c9be9SBarry Smith
12297f5c9be9SBarry Smith Input Parameter:
12302fe279fdSBarry Smith . tao - the `Tao` context
12317f5c9be9SBarry Smith
12327f5c9be9SBarry Smith Output Parameter:
12332fe279fdSBarry Smith . admm_tao - the parent `Tao` context
12346285c0a3SHansol Suh
12356285c0a3SHansol Suh Level: advanced
1236a5a2f7acSBarry Smith
1237db781477SPatrick Sanan .seealso: `TAOADMM`
12386285c0a3SHansol Suh @*/
TaoGetADMMParentTao(Tao tao,Tao * admm_tao)1239d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao)
1240d71ae5a4SJacob Faibussowitsch {
12416285c0a3SHansol Suh PetscFunctionBegin;
12426285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12439566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)tao, "TaoGetADMMParentTao_ADMM", (PetscObject *)admm_tao));
12443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12456285c0a3SHansol Suh }
12466285c0a3SHansol Suh
12476285c0a3SHansol Suh /*@
12482fe279fdSBarry Smith TaoADMMGetDualVector - Returns the dual vector associated with the current `TAOADMM` state
12496285c0a3SHansol Suh
12506285c0a3SHansol Suh Not Collective
12516285c0a3SHansol Suh
12526285c0a3SHansol Suh Input Parameter:
12532fe279fdSBarry Smith . tao - the `Tao` context
12547f5c9be9SBarry Smith
12557f5c9be9SBarry Smith Output Parameter:
12567f5c9be9SBarry Smith . Y - the current solution
12576285c0a3SHansol Suh
12586285c0a3SHansol Suh Level: intermediate
1259a5a2f7acSBarry Smith
1260db781477SPatrick Sanan .seealso: `TAOADMM`
12616285c0a3SHansol Suh @*/
TaoADMMGetDualVector(Tao tao,Vec * Y)1262d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y)
1263d71ae5a4SJacob Faibussowitsch {
12646285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM *)tao->data;
12656285c0a3SHansol Suh
12666285c0a3SHansol Suh PetscFunctionBegin;
12676285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12686285c0a3SHansol Suh *Y = am->y;
12693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12706285c0a3SHansol Suh }
12716285c0a3SHansol Suh
12726285c0a3SHansol Suh /*@
12732fe279fdSBarry Smith TaoADMMSetRegularizerType - Set regularizer type for `TAOADMM` routine
12746285c0a3SHansol Suh
12756285c0a3SHansol Suh Not Collective
12766285c0a3SHansol Suh
1277d8d19677SJose E. Roman Input Parameters:
12782fe279fdSBarry Smith + tao - the `Tao` context
12797f5c9be9SBarry Smith - type - regularizer type
12807f5c9be9SBarry Smith
12813c7db156SBarry Smith Options Database Key:
1282147403d9SBarry Smith . -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer
12836285c0a3SHansol Suh
12846285c0a3SHansol Suh Level: intermediate
12856285c0a3SHansol Suh
1286db781477SPatrick Sanan .seealso: `TaoADMMGetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
12876285c0a3SHansol Suh @*/
TaoADMMSetRegularizerType(Tao tao,TaoADMMRegularizerType type)1288d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type)
1289d71ae5a4SJacob Faibussowitsch {
12906285c0a3SHansol Suh PetscFunctionBegin;
12916285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12926285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao, type, 2);
1293cac4c232SBarry Smith PetscTryMethod(tao, "TaoADMMSetRegularizerType_C", (Tao, TaoADMMRegularizerType), (tao, type));
12943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12956285c0a3SHansol Suh }
12966285c0a3SHansol Suh
12976285c0a3SHansol Suh /*@
12982fe279fdSBarry Smith TaoADMMGetRegularizerType - Gets the type of regularizer routine for `TAOADMM`
12996285c0a3SHansol Suh
13006285c0a3SHansol Suh Not Collective
13016285c0a3SHansol Suh
13026285c0a3SHansol Suh Input Parameter:
13032fe279fdSBarry Smith . tao - the `Tao` context
13046285c0a3SHansol Suh
13056285c0a3SHansol Suh Output Parameter:
13066285c0a3SHansol Suh . type - the type of regularizer
13076285c0a3SHansol Suh
13086285c0a3SHansol Suh Level: intermediate
13096285c0a3SHansol Suh
1310db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
13116285c0a3SHansol Suh @*/
TaoADMMGetRegularizerType(Tao tao,TaoADMMRegularizerType * type)1312d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type)
1313d71ae5a4SJacob Faibussowitsch {
13146285c0a3SHansol Suh PetscFunctionBegin;
13156285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1316cac4c232SBarry Smith PetscUseMethod(tao, "TaoADMMGetRegularizerType_C", (Tao, TaoADMMRegularizerType *), (tao, type));
13173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13186285c0a3SHansol Suh }
13196285c0a3SHansol Suh
13206285c0a3SHansol Suh /*@
13212fe279fdSBarry Smith TaoADMMSetUpdateType - Set update routine for `TAOADMM` routine
13226285c0a3SHansol Suh
13236285c0a3SHansol Suh Not Collective
13246285c0a3SHansol Suh
1325d8d19677SJose E. Roman Input Parameters:
13262fe279fdSBarry Smith + tao - the `Tao` context
13276285c0a3SHansol Suh - type - spectral parameter update type
13286285c0a3SHansol Suh
13296285c0a3SHansol Suh Level: intermediate
13306285c0a3SHansol Suh
1331db781477SPatrick Sanan .seealso: `TaoADMMGetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
13326285c0a3SHansol Suh @*/
TaoADMMSetUpdateType(Tao tao,TaoADMMUpdateType type)1333d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type)
1334d71ae5a4SJacob Faibussowitsch {
13356285c0a3SHansol Suh PetscFunctionBegin;
13366285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
13376285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao, type, 2);
1338cac4c232SBarry Smith PetscTryMethod(tao, "TaoADMMSetUpdateType_C", (Tao, TaoADMMUpdateType), (tao, type));
13393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13406285c0a3SHansol Suh }
13416285c0a3SHansol Suh
13426285c0a3SHansol Suh /*@
13432fe279fdSBarry Smith TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for `TAOADMM`
13446285c0a3SHansol Suh
13456285c0a3SHansol Suh Not Collective
13466285c0a3SHansol Suh
13476285c0a3SHansol Suh Input Parameter:
13482fe279fdSBarry Smith . tao - the `Tao` context
13496285c0a3SHansol Suh
13506285c0a3SHansol Suh Output Parameter:
13516285c0a3SHansol Suh . type - the type of spectral penalty update routine
13526285c0a3SHansol Suh
13536285c0a3SHansol Suh Level: intermediate
13546285c0a3SHansol Suh
1355db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
13566285c0a3SHansol Suh @*/
TaoADMMGetUpdateType(Tao tao,TaoADMMUpdateType * type)1357d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type)
1358d71ae5a4SJacob Faibussowitsch {
13596285c0a3SHansol Suh PetscFunctionBegin;
13606285c0a3SHansol Suh PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1361cac4c232SBarry Smith PetscUseMethod(tao, "TaoADMMGetUpdateType_C", (Tao, TaoADMMUpdateType *), (tao, type));
13623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13636285c0a3SHansol Suh }
1364