xref: /petsc/src/tao/constrained/impls/admm/admm.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 #include <../src/tao/constrained/impls/admm/admm.h> /*I "petsctao.h" I*/
2 #include <petsctao.h>
3 #include <petsc/private/petscimpl.h>
4 
5 /* Updates terminating criteria
6  *
7  * 1  ||r_k|| = ||Ax+Bz-c|| =< catol_admm* max{||Ax||,||Bz||,||c||}
8  *
9  * 2. Updates dual residual, d_k
10  *
11  * 3. ||d_k|| = ||mu*A^T*B(z_k-z_{k-1})|| =< gatol_admm * ||A^Ty||   */
12 
13 static PetscBool cited = PETSC_FALSE;
14 static const char citation[] =
15   "@misc{xu2017adaptive,\n"
16   "   title={Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation},\n"
17   "   author={Zheng Xu and Mario A. T. Figueiredo and Xiaoming Yuan and Christoph Studer and Tom Goldstein},\n"
18   "   year={2017},\n"
19   "   eprint={1704.02712},\n"
20   "   archivePrefix={arXiv},\n"
21   "   primaryClass={cs.CV}\n"
22   "}  \n";
23 
24 const char *const TaoADMMRegularizerTypes[] = {"REGULARIZER_USER","REGULARIZER_SOFT_THRESH","TaoADMMRegularizerType","TAO_ADMM_",NULL};
25 const char *const TaoADMMUpdateTypes[]      = {"UPDATE_BASIC","UPDATE_ADAPTIVE","UPDATE_ADAPTIVE_RELAXED","TaoADMMUpdateType","TAO_ADMM_",NULL};
26 const char *const TaoALMMTypes[]            = {"CLASSIC","PHR","TaoALMMType","TAO_ALMM_",NULL};
27 
28 static PetscErrorCode TaoADMMToleranceUpdate(Tao tao)
29 {
30   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
31   PetscReal      Axnorm,Bznorm,ATynorm,temp;
32   Vec            tempJR,tempL;
33   Tao            mis;
34 
35   PetscFunctionBegin;
36   mis    = am->subsolverX;
37   tempJR = am->workJacobianRight;
38   tempL  = am->workLeft;
39   /* ATy */
40   PetscCall(TaoComputeJacobianEquality(mis, am->y, mis->jacobian_equality, mis->jacobian_equality_pre));
41   PetscCall(MatMultTranspose(mis->jacobian_equality,am->y,tempJR));
42   PetscCall(VecNorm(tempJR,NORM_2,&ATynorm));
43   /* dualres = mu * ||AT(Bz-Bzold)||_2 */
44   PetscCall(VecWAXPY(tempJR,-1.,am->Bzold,am->Bz));
45   PetscCall(MatMultTranspose(mis->jacobian_equality,tempJR,tempL));
46   PetscCall(VecNorm(tempL,NORM_2,&am->dualres));
47   am->dualres *= am->mu;
48 
49   /* ||Ax||_2, ||Bz||_2 */
50   PetscCall(VecNorm(am->Ax,NORM_2,&Axnorm));
51   PetscCall(VecNorm(am->Bz,NORM_2,&Bznorm));
52 
53   /* Set catol to be catol_admm *  max{||Ax||,||Bz||,||c||} *
54    * Set gatol to be gatol_admm *  ||A^Ty|| *
55    * while cnorm is ||r_k||_2, and gnorm is ||d_k||_2 */
56   temp = am->catol_admm * PetscMax(Axnorm, (!am->const_norm) ? Bznorm : PetscMax(Bznorm,am->const_norm));
57   PetscCall(TaoSetConstraintTolerances(tao,temp,PETSC_DEFAULT));
58   PetscCall(TaoSetTolerances(tao, am->gatol_admm*ATynorm, PETSC_DEFAULT,PETSC_DEFAULT));
59   PetscFunctionReturn(0);
60 }
61 
62 /* Penaly Update for Adaptive ADMM. */
63 static PetscErrorCode AdaptiveADMMPenaltyUpdate(Tao tao)
64 {
65   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
66   PetscReal      ydiff_norm, yhatdiff_norm, Axdiff_norm, Bzdiff_norm, Axyhat, Bzy, a_sd, a_mg, a_k, b_sd, b_mg, b_k;
67   PetscBool      hflag, gflag;
68   Vec            tempJR,tempJR2;
69 
70   PetscFunctionBegin;
71   tempJR  = am->workJacobianRight;
72   tempJR2 = am->workJacobianRight2;
73   hflag   = PETSC_FALSE;
74   gflag   = PETSC_FALSE;
75   a_k     = -1;
76   b_k     = -1;
77 
78   PetscCall(VecWAXPY(tempJR,-1.,am->Axold,am->Ax));
79   PetscCall(VecWAXPY(tempJR2,-1.,am->yhatold,am->yhat));
80   PetscCall(VecNorm(tempJR,NORM_2,&Axdiff_norm));
81   PetscCall(VecNorm(tempJR2,NORM_2,&yhatdiff_norm));
82   PetscCall(VecDot(tempJR,tempJR2,&Axyhat));
83 
84   PetscCall(VecWAXPY(tempJR,-1.,am->Bz0,am->Bz));
85   PetscCall(VecWAXPY(tempJR2,-1.,am->y,am->y0));
86   PetscCall(VecNorm(tempJR,NORM_2,&Bzdiff_norm));
87   PetscCall(VecNorm(tempJR2,NORM_2,&ydiff_norm));
88   PetscCall(VecDot(tempJR,tempJR2,&Bzy));
89 
90   if (Axyhat > am->orthval*Axdiff_norm*yhatdiff_norm + am->mueps) {
91     hflag = PETSC_TRUE;
92     a_sd  = PetscSqr(yhatdiff_norm)/Axyhat; /* alphaSD */
93     a_mg  = Axyhat/PetscSqr(Axdiff_norm);   /* alphaMG */
94     a_k   = (a_mg/a_sd) > 0.5 ? a_mg : a_sd - 0.5*a_mg;
95   }
96   if (Bzy > am->orthval*Bzdiff_norm*ydiff_norm + am->mueps) {
97     gflag = PETSC_TRUE;
98     b_sd  = PetscSqr(ydiff_norm)/Bzy;  /* betaSD */
99     b_mg  = Bzy/PetscSqr(Bzdiff_norm); /* betaMG */
100     b_k   = (b_mg/b_sd) > 0.5 ? b_mg : b_sd - 0.5*b_mg;
101   }
102   am->muold = am->mu;
103   if (gflag && hflag) {
104     am->mu = PetscSqrtReal(a_k*b_k);
105   } else if (hflag) {
106     am->mu = a_k;
107   } else if (gflag) {
108     am->mu = b_k;
109   }
110   if (am->mu > am->muold) {
111     am->mu = am->muold;
112   }
113   if (am->mu < am->mumin) {
114     am->mu = am->mumin;
115   }
116   PetscFunctionReturn(0);
117 }
118 
119 static PetscErrorCode  TaoADMMSetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType type)
120 {
121   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
122 
123   PetscFunctionBegin;
124   am->regswitch = type;
125   PetscFunctionReturn(0);
126 }
127 
128 static PetscErrorCode  TaoADMMGetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType *type)
129 {
130   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
131 
132   PetscFunctionBegin;
133   *type = am->regswitch;
134   PetscFunctionReturn(0);
135 }
136 
137 static PetscErrorCode  TaoADMMSetUpdateType_ADMM(Tao tao, TaoADMMUpdateType type)
138 {
139   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
140 
141   PetscFunctionBegin;
142   am->update = type;
143   PetscFunctionReturn(0);
144 }
145 
146 static PetscErrorCode  TaoADMMGetUpdateType_ADMM(Tao tao, TaoADMMUpdateType *type)
147 {
148   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
149 
150   PetscFunctionBegin;
151   *type = am->update;
152   PetscFunctionReturn(0);
153 }
154 
155 /* This routine updates Jacobians with new x,z vectors,
156  * and then updates Ax and Bz vectors, then computes updated residual vector*/
157 static PetscErrorCode ADMMUpdateConstraintResidualVector(Tao tao, Vec x, Vec z, Vec Ax, Vec Bz, Vec residual)
158 {
159   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
160   Tao            mis,reg;
161 
162   PetscFunctionBegin;
163   mis  = am->subsolverX;
164   reg  = am->subsolverZ;
165   PetscCall(TaoComputeJacobianEquality(mis, x, mis->jacobian_equality, mis->jacobian_equality_pre));
166   PetscCall(MatMult(mis->jacobian_equality,x,Ax));
167   PetscCall(TaoComputeJacobianEquality(reg, z, reg->jacobian_equality, reg->jacobian_equality_pre));
168   PetscCall(MatMult(reg->jacobian_equality,z,Bz));
169 
170   PetscCall(VecWAXPY(residual,1.,Bz,Ax));
171   if (am->constraint != NULL) {
172     PetscCall(VecAXPY(residual,-1.,am->constraint));
173   }
174   PetscFunctionReturn(0);
175 }
176 
177 /* Updates Augmented Lagrangians to given routines *
178  * For subsolverX, routine needs to be ComputeObjectiveAndGraidnet
179  * Separate Objective and Gradient routines are not supported.  */
180 static PetscErrorCode SubObjGradUpdate(Tao tao, Vec x, PetscReal *f, Vec g, void *ptr)
181 {
182   Tao            parent = (Tao)ptr;
183   TAO_ADMM       *am    = (TAO_ADMM*)parent->data;
184   PetscReal      temp,temp2;
185   Vec            tempJR;
186 
187   PetscFunctionBegin;
188   tempJR = am->workJacobianRight;
189   PetscCall(ADMMUpdateConstraintResidualVector(parent, x, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
190   PetscCall((*am->ops->misfitobjgrad)(am->subsolverX,x,f,g,am->misfitobjgradP));
191 
192   am->last_misfit_val = *f;
193   /* Objective  Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */
194   PetscCall(VecTDot(am->residual,am->y,&temp));
195   PetscCall(VecTDot(am->residual,am->residual,&temp2));
196   *f   += temp + (am->mu/2)*temp2;
197 
198   /* Gradient. Add + mu*AT(Ax+Bz-c) + yTA*/
199   PetscCall(MatMultTranspose(tao->jacobian_equality,am->residual,tempJR));
200   PetscCall(VecAXPY(g,am->mu,tempJR));
201   PetscCall(MatMultTranspose(tao->jacobian_equality,am->y,tempJR));
202   PetscCall(VecAXPY(g,1.,tempJR));
203   PetscFunctionReturn(0);
204 }
205 
206 /* Updates Augmented Lagrangians to given routines
207  * For subsolverZ, routine needs to be ComputeObjectiveAndGraidnet
208  * Separate Objective and Gradient routines are not supported.  */
209 static PetscErrorCode RegObjGradUpdate(Tao tao, Vec z, PetscReal *f, Vec g, void *ptr)
210 {
211   Tao            parent = (Tao)ptr;
212   TAO_ADMM       *am    = (TAO_ADMM*)parent->data;
213   PetscReal      temp,temp2;
214   Vec            tempJR;
215 
216   PetscFunctionBegin;
217   tempJR = am->workJacobianRight;
218   PetscCall(ADMMUpdateConstraintResidualVector(parent, am->subsolverX->solution, z, am->Ax, am->Bz, am->residual));
219   PetscCall((*am->ops->regobjgrad)(am->subsolverZ,z,f,g,am->regobjgradP));
220   am->last_reg_val= *f;
221   /* Objective  Add  + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */
222   PetscCall(VecTDot(am->residual,am->y,&temp));
223   PetscCall(VecTDot(am->residual,am->residual,&temp2));
224   *f   += temp + (am->mu/2)*temp2;
225 
226   /* Gradient. Add + mu*BT(Ax+Bz-c) + yTB*/
227   PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality,am->residual,tempJR));
228   PetscCall(VecAXPY(g,am->mu,tempJR));
229   PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality,am->y,tempJR));
230   PetscCall(VecAXPY(g,1.,tempJR));
231   PetscFunctionReturn(0);
232 }
233 
234 /* Computes epsilon padded L1 norm lambda*sum(sqrt(x^2+eps^2)-eps */
235 static PetscErrorCode ADMML1EpsilonNorm(Tao tao, Vec x, PetscReal eps, PetscReal *norm)
236 {
237   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
238   PetscInt       N;
239 
240   PetscFunctionBegin;
241   PetscCall(VecGetSize(am->workLeft,&N));
242   PetscCall(VecPointwiseMult(am->workLeft,x,x));
243   PetscCall(VecShift(am->workLeft,am->l1epsilon*am->l1epsilon));
244   PetscCall(VecSqrtAbs(am->workLeft));
245   PetscCall(VecSum(am->workLeft,norm));
246   *norm += N*am->l1epsilon;
247   *norm *= am->lambda;
248   PetscFunctionReturn(0);
249 }
250 
251 static PetscErrorCode ADMMInternalHessianUpdate(Mat H, Mat Constraint, PetscBool Identity, void *ptr)
252 {
253   TAO_ADMM       *am = (TAO_ADMM*)ptr;
254 
255   PetscFunctionBegin;
256   switch (am->update) {
257   case (TAO_ADMM_UPDATE_BASIC):
258     break;
259   case (TAO_ADMM_UPDATE_ADAPTIVE):
260   case (TAO_ADMM_UPDATE_ADAPTIVE_RELAXED):
261     if (H && (am->muold != am->mu)) {
262       if (!Identity) {
263         PetscCall(MatAXPY(H,am->mu-am->muold,Constraint,DIFFERENT_NONZERO_PATTERN));
264       } else {
265         PetscCall(MatShift(H,am->mu-am->muold));
266       }
267     }
268     break;
269   }
270   PetscFunctionReturn(0);
271 }
272 
273 /* Updates Hessian - adds second derivative of augmented Lagrangian
274  * H \gets H + \rho*ATA
275  * Here, \rho does not change in TAO_ADMM_UPDATE_BASIC - thus no-op
276  * For ADAPTAIVE,ADAPTIVE_RELAXED,
277  * H \gets H + (\rho-\rhoold)*ATA
278  * Here, we assume that A is linear constraint i.e., doesnt change.
279  * Thus, for both ADAPTIVE, and RELAXED, ATA matrix is pre-set (except for A=I (null case)) see TaoSetUp_ADMM */
280 static PetscErrorCode SubHessianUpdate(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
281 {
282   Tao            parent = (Tao)ptr;
283   TAO_ADMM       *am    = (TAO_ADMM*)parent->data;
284 
285   PetscFunctionBegin;
286   if (am->Hxchange) {
287     /* Case where Hessian gets updated with respect to x vector input. */
288     PetscCall((*am->ops->misfithess)(am->subsolverX,x,H,Hpre,am->misfithessP));
289     PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am));
290   } else if (am->Hxbool) {
291     /* Hessian doesn't get updated. H(x) = c */
292     /* Update Lagrangian only only per TAO call */
293     PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am));
294     am->Hxbool = PETSC_FALSE;
295   }
296   PetscFunctionReturn(0);
297 }
298 
299 /* Same as SubHessianUpdate, except for B matrix instead of A matrix */
300 static PetscErrorCode RegHessianUpdate(Tao tao, Vec z, Mat H, Mat Hpre, void *ptr)
301 {
302   Tao            parent = (Tao)ptr;
303   TAO_ADMM       *am    = (TAO_ADMM*)parent->data;
304 
305   PetscFunctionBegin;
306 
307   if (am->Hzchange) {
308     /* Case where Hessian gets updated with respect to x vector input. */
309     PetscCall((*am->ops->reghess)(am->subsolverZ,z,H,Hpre,am->reghessP));
310     PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am));
311   } else if (am->Hzbool) {
312     /* Hessian doesn't get updated. H(x) = c */
313     /* Update Lagrangian only only per TAO call */
314     PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am));
315     am->Hzbool = PETSC_FALSE;
316   }
317   PetscFunctionReturn(0);
318 }
319 
320 /* Shell Matrix routine for A matrix.
321  * This gets used when user puts NULL for
322  * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...)
323  * Essentially sets A=I*/
324 static PetscErrorCode JacobianIdentity(Mat mat,Vec in,Vec out)
325 {
326   PetscFunctionBegin;
327   PetscCall(VecCopy(in,out));
328   PetscFunctionReturn(0);
329 }
330 
331 /* Shell Matrix routine for B matrix.
332  * This gets used when user puts NULL for
333  * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...)
334  * Sets B=-I */
335 static PetscErrorCode JacobianIdentityB(Mat mat,Vec in,Vec out)
336 {
337   PetscFunctionBegin;
338   PetscCall(VecCopy(in,out));
339   PetscCall(VecScale(out,-1.));
340   PetscFunctionReturn(0);
341 }
342 
343 /* Solve f(x) + g(z) s.t. Ax + Bz = c */
344 static PetscErrorCode TaoSolve_ADMM(Tao tao)
345 {
346   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
347   PetscInt       N;
348   PetscReal      reg_func;
349   PetscBool      is_reg_shell;
350   Vec            tempL;
351 
352   PetscFunctionBegin;
353   if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) {
354     PetscCheck(am->subsolverX->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetMisfitConstraintJacobian() first");
355     PetscCheck(am->subsolverZ->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerConstraintJacobian() first");
356     if (am->constraint != NULL) {
357       PetscCall(VecNorm(am->constraint,NORM_2,&am->const_norm));
358     }
359   }
360   tempL = am->workLeft;
361   PetscCall(VecGetSize(tempL,&N));
362 
363   if (am->Hx && am->ops->misfithess) {
364     PetscCall(TaoSetHessian(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao));
365   }
366 
367   if (!am->zJI) {
368     /* Currently, B is assumed to be a linear system, i.e., not getting updated*/
369     PetscCall(MatTransposeMatMult(am->JB,am->JB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->BTB)));
370   }
371   if (!am->xJI) {
372     /* Currently, A is assumed to be a linear system, i.e., not getting updated*/
373     PetscCall(MatTransposeMatMult(am->subsolverX->jacobian_equality,am->subsolverX->jacobian_equality,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->ATA)));
374   }
375 
376   is_reg_shell = PETSC_FALSE;
377 
378   PetscCall(PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell));
379 
380   if (!is_reg_shell) {
381     switch (am->regswitch) {
382     case (TAO_ADMM_REGULARIZER_USER):
383       PetscCheck(am->ops->regobjgrad,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerObjectiveAndGradientRoutine() first if one wishes to use TAO_ADMM_REGULARIZER_USER with non-TAOSHELL type");
384       break;
385     case (TAO_ADMM_REGULARIZER_SOFT_THRESH):
386       /* Soft Threshold. */
387       break;
388     }
389     if (am->ops->regobjgrad) PetscCall(TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao));
390     if (am->Hz && am->ops->reghess) {
391       PetscCall(TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao));
392     }
393   }
394 
395   switch (am->update) {
396   case TAO_ADMM_UPDATE_BASIC:
397     if (am->subsolverX->hessian) {
398       /* In basic case, Hessian does not get updated w.r.t. to spectral penalty
399        * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/
400       if (!am->xJI) {
401         PetscCall(MatAXPY(am->subsolverX->hessian,am->mu,am->ATA,DIFFERENT_NONZERO_PATTERN));
402       } else {
403         PetscCall(MatShift(am->subsolverX->hessian,am->mu));
404       }
405     }
406     if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) {
407       if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) {
408         PetscCall(MatAXPY(am->subsolverZ->hessian,am->mu,am->BTB,DIFFERENT_NONZERO_PATTERN));
409       } else {
410         PetscCall(MatShift(am->subsolverZ->hessian,am->mu));
411       }
412     }
413     break;
414   case TAO_ADMM_UPDATE_ADAPTIVE:
415   case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
416     break;
417   }
418 
419   PetscCall(PetscCitationsRegister(citation,&cited));
420   tao->reason = TAO_CONTINUE_ITERATING;
421 
422   while (tao->reason == TAO_CONTINUE_ITERATING) {
423     if (tao->ops->update) PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
424     PetscCall(VecCopy(am->Bz, am->Bzold));
425 
426     /* x update */
427     PetscCall(TaoSolve(am->subsolverX));
428     PetscCall(TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre));
429     PetscCall(MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution,am->Ax));
430 
431     am->Hxbool = PETSC_TRUE;
432 
433     /* z update */
434     switch (am->regswitch) {
435     case TAO_ADMM_REGULARIZER_USER:
436       PetscCall(TaoSolve(am->subsolverZ));
437       break;
438     case TAO_ADMM_REGULARIZER_SOFT_THRESH:
439       /* L1 assumes A,B jacobians are identity nxn matrix */
440       PetscCall(VecWAXPY(am->workJacobianRight,1/am->mu,am->y,am->Ax));
441       PetscCall(TaoSoftThreshold(am->workJacobianRight,-am->lambda/am->mu,am->lambda/am->mu,am->subsolverZ->solution));
442       break;
443     }
444     am->Hzbool = PETSC_TRUE;
445     /* Returns Ax + Bz - c with updated Ax,Bz vectors */
446     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
447     /* Dual variable, y += y + mu*(Ax+Bz-c) */
448     PetscCall(VecWAXPY(am->y, am->mu, am->residual, am->yold));
449 
450     /* stopping tolerance update */
451     PetscCall(TaoADMMToleranceUpdate(tao));
452 
453     /* Updating Spectral Penalty */
454     switch (am->update) {
455     case TAO_ADMM_UPDATE_BASIC:
456       am->muold = am->mu;
457       break;
458     case TAO_ADMM_UPDATE_ADAPTIVE:
459     case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
460       if (tao->niter == 0) {
461         PetscCall(VecCopy(am->y, am->y0));
462         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
463         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
464         PetscCall(VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold));
465         PetscCall(VecCopy(am->Ax, am->Axold));
466         PetscCall(VecCopy(am->Bz, am->Bz0));
467         am->muold = am->mu;
468       } else if (tao->niter % am->T == 1) {
469         /* we have compute Bzold in a previous iteration, and we computed Ax above */
470         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
471         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
472         PetscCall(VecWAXPY(am->yhat, -am->mu, am->residual, am->yold));
473         PetscCall(AdaptiveADMMPenaltyUpdate(tao));
474         PetscCall(VecCopy(am->Ax, am->Axold));
475         PetscCall(VecCopy(am->Bz, am->Bz0));
476         PetscCall(VecCopy(am->yhat, am->yhatold));
477         PetscCall(VecCopy(am->y, am->y0));
478       } else {
479         am->muold = am->mu;
480       }
481       break;
482     default:
483       break;
484     }
485     tao->niter++;
486 
487     /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/
488     switch (am->regswitch) {
489     case TAO_ADMM_REGULARIZER_USER:
490       if (is_reg_shell) {
491         PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,&reg_func));
492       } else {
493         (*am->ops->regobjgrad)(am->subsolverZ,am->subsolverX->solution,&reg_func,tempL,am->regobjgradP);
494       }
495       break;
496     case TAO_ADMM_REGULARIZER_SOFT_THRESH:
497       PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,&reg_func));
498       break;
499     }
500     PetscCall(VecCopy(am->y,am->yold));
501     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
502     PetscCall(VecNorm(am->residual,NORM_2,&am->resnorm));
503     PetscCall(TaoLogConvergenceHistory(tao,am->last_misfit_val + reg_func,am->dualres,am->resnorm,tao->ksp_its));
504 
505     PetscCall(TaoMonitor(tao,tao->niter,am->last_misfit_val + reg_func,am->dualres,am->resnorm,1.0));
506     PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
507   }
508   /* Update vectors */
509   PetscCall(VecCopy(am->subsolverX->solution,tao->solution));
510   PetscCall(VecCopy(am->subsolverX->gradient,tao->gradient));
511   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", NULL));
512   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", NULL));
513   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL));
514   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL));
515   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL));
516   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL));
517   PetscFunctionReturn(0);
518 }
519 
520 static PetscErrorCode TaoSetFromOptions_ADMM(PetscOptionItems *PetscOptionsObject,Tao tao)
521 {
522   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
523 
524   PetscFunctionBegin;
525   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. ");
526   PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient","regularizer constant","",am->lambda,&am->lambda,NULL));
527   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty","Constant for Augmented Lagrangian term.","",am->mu,&am->mu,NULL));
528   PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter","x relaxation parameter for Z update.","",am->gamma,&am->gamma,NULL));
529   PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor","ADMM dynamic tolerance update factor.","",am->tol,&am->tol,NULL));
530   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor","ADMM spectral penalty update curvature safeguard value.","",am->orthval,&am->orthval,NULL));
531   PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty","Set ADMM minimum spectral penalty.","",am->mumin,&am->mumin,NULL));
532   PetscCall(PetscOptionsEnum("-tao_admm_dual_update","Lagrangian dual update policy","TaoADMMUpdateType",TaoADMMUpdateTypes,(PetscEnum)am->update,(PetscEnum*)&am->update,NULL));
533   PetscCall(PetscOptionsEnum("-tao_admm_regularizer_type","ADMM regularizer update rule","TaoADMMRegularizerType",TaoADMMRegularizerTypes,(PetscEnum)am->regswitch,(PetscEnum*)&am->regswitch,NULL));
534   PetscOptionsHeadEnd();
535   PetscCall(TaoSetFromOptions(am->subsolverX));
536   if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) {
537     PetscCall(TaoSetFromOptions(am->subsolverZ));
538   }
539   PetscFunctionReturn(0);
540 }
541 
542 static PetscErrorCode TaoView_ADMM(Tao tao,PetscViewer viewer)
543 {
544   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
545 
546   PetscFunctionBegin;
547   PetscCall(PetscViewerASCIIPushTab(viewer));
548   PetscCall(TaoView(am->subsolverX,viewer));
549   PetscCall(TaoView(am->subsolverZ,viewer));
550   PetscCall(PetscViewerASCIIPopTab(viewer));
551   PetscFunctionReturn(0);
552 }
553 
554 static PetscErrorCode TaoSetUp_ADMM(Tao tao)
555 {
556   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
557   PetscInt       n,N,M;
558 
559   PetscFunctionBegin;
560   PetscCall(VecGetLocalSize(tao->solution,&n));
561   PetscCall(VecGetSize(tao->solution,&N));
562   /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */
563   if (!am->JB) {
564     am->zJI   = PETSC_TRUE;
565     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JB));
566     PetscCall(MatShellSetOperation(am->JB,MATOP_MULT,(void (*)(void))JacobianIdentityB));
567     PetscCall(MatShellSetOperation(am->JB,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentityB));
568     am->JBpre = am->JB;
569   }
570   if (!am->JA) {
571     am->xJI   = PETSC_TRUE;
572     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JA));
573     PetscCall(MatShellSetOperation(am->JA,MATOP_MULT,(void (*)(void))JacobianIdentity));
574     PetscCall(MatShellSetOperation(am->JA,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentity));
575     am->JApre = am->JA;
576   }
577   PetscCall(MatCreateVecs(am->JA,NULL,&am->Ax));
578   if (!tao->gradient) {
579     PetscCall(VecDuplicate(tao->solution,&tao->gradient));
580   }
581   PetscCall(TaoSetSolution(am->subsolverX, tao->solution));
582   if (!am->z) {
583     PetscCall(VecDuplicate(tao->solution,&am->z));
584     PetscCall(VecSet(am->z,0.0));
585   }
586   PetscCall(TaoSetSolution(am->subsolverZ, am->z));
587   if (!am->workLeft) {
588     PetscCall(VecDuplicate(tao->solution,&am->workLeft));
589   }
590   if (!am->Axold) {
591     PetscCall(VecDuplicate(am->Ax,&am->Axold));
592   }
593   if (!am->workJacobianRight) {
594     PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight));
595   }
596   if (!am->workJacobianRight2) {
597     PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight2));
598   }
599   if (!am->Bz) {
600     PetscCall(VecDuplicate(am->Ax,&am->Bz));
601   }
602   if (!am->Bzold) {
603     PetscCall(VecDuplicate(am->Ax,&am->Bzold));
604   }
605   if (!am->Bz0) {
606     PetscCall(VecDuplicate(am->Ax,&am->Bz0));
607   }
608   if (!am->y) {
609     PetscCall(VecDuplicate(am->Ax,&am->y));
610     PetscCall(VecSet(am->y,0.0));
611   }
612   if (!am->yold) {
613     PetscCall(VecDuplicate(am->Ax,&am->yold));
614     PetscCall(VecSet(am->yold,0.0));
615   }
616   if (!am->y0) {
617     PetscCall(VecDuplicate(am->Ax,&am->y0));
618     PetscCall(VecSet(am->y0,0.0));
619   }
620   if (!am->yhat) {
621     PetscCall(VecDuplicate(am->Ax,&am->yhat));
622     PetscCall(VecSet(am->yhat,0.0));
623   }
624   if (!am->yhatold) {
625     PetscCall(VecDuplicate(am->Ax,&am->yhatold));
626     PetscCall(VecSet(am->yhatold,0.0));
627   }
628   if (!am->residual) {
629     PetscCall(VecDuplicate(am->Ax,&am->residual));
630     PetscCall(VecSet(am->residual,0.0));
631   }
632   if (!am->constraint) {
633     am->constraint = NULL;
634   } else {
635     PetscCall(VecGetSize(am->constraint,&M));
636     PetscCheck(M == N,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Solution vector and constraint vector must be of same size!");
637   }
638 
639   /* Save changed tao tolerance for adaptive tolerance */
640   if (tao->gatol_changed) {
641     am->gatol_admm = tao->gatol;
642   }
643   if (tao->catol_changed) {
644     am->catol_admm = tao->catol;
645   }
646 
647   /*Update spectral and dual elements to X subsolver */
648   PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao));
649   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX,am->JA,am->JApre, am->ops->misfitjac, am->misfitjacobianP));
650   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ,am->JB,am->JBpre, am->ops->regjac, am->regjacobianP));
651   PetscFunctionReturn(0);
652 }
653 
654 static PetscErrorCode TaoDestroy_ADMM(Tao tao)
655 {
656   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
657 
658   PetscFunctionBegin;
659   PetscCall(VecDestroy(&am->z));
660   PetscCall(VecDestroy(&am->Ax));
661   PetscCall(VecDestroy(&am->Axold));
662   PetscCall(VecDestroy(&am->Bz));
663   PetscCall(VecDestroy(&am->Bzold));
664   PetscCall(VecDestroy(&am->Bz0));
665   PetscCall(VecDestroy(&am->residual));
666   PetscCall(VecDestroy(&am->y));
667   PetscCall(VecDestroy(&am->yold));
668   PetscCall(VecDestroy(&am->y0));
669   PetscCall(VecDestroy(&am->yhat));
670   PetscCall(VecDestroy(&am->yhatold));
671   PetscCall(VecDestroy(&am->workLeft));
672   PetscCall(VecDestroy(&am->workJacobianRight));
673   PetscCall(VecDestroy(&am->workJacobianRight2));
674 
675   PetscCall(MatDestroy(&am->JA));
676   PetscCall(MatDestroy(&am->JB));
677   if (!am->xJI) {
678     PetscCall(MatDestroy(&am->JApre));
679   }
680   if (!am->zJI) {
681     PetscCall(MatDestroy(&am->JBpre));
682   }
683   if (am->Hx) {
684     PetscCall(MatDestroy(&am->Hx));
685     PetscCall(MatDestroy(&am->Hxpre));
686   }
687   if (am->Hz) {
688     PetscCall(MatDestroy(&am->Hz));
689     PetscCall(MatDestroy(&am->Hzpre));
690   }
691   PetscCall(MatDestroy(&am->ATA));
692   PetscCall(MatDestroy(&am->BTB));
693   PetscCall(TaoDestroy(&am->subsolverX));
694   PetscCall(TaoDestroy(&am->subsolverZ));
695   am->parent = NULL;
696   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL));
697   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL));
698   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL));
699   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL));
700   PetscCall(PetscFree(tao->data));
701   PetscFunctionReturn(0);
702 }
703 
704 /*MC
705 
706   TAOADMM - Alternating direction method of multipliers method fo solving linear problems with
707             constraints. in a min_x f(x) + g(z)  s.t. Ax+Bz=c.
708             This algorithm employs two sub Tao solvers, of which type can be specified
709             by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers.
710             Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via
711             TaoADMMSet{Misfit,Regularizer}HessianChangeStatus.
712             Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type.
713             There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update,
714             currently there are basic option and adaptive option.
715             Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian.
716             c can be set with TaoADMMSetConstraintVectorRHS.
717             The user can also provide regularizer weight for second subsolver.
718 
719   References:
720 . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom
721           "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation"
722           The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017.
723 
724   Options Database Keys:
725 + -tao_admm_regularizer_coefficient        - regularizer constant (default 1.e-6)
726 . -tao_admm_spectral_penalty               - Constant for Augmented Lagrangian term (default 1.)
727 . -tao_admm_relaxation_parameter           - relaxation parameter for Z update (default 1.)
728 . -tao_admm_tolerance_update_factor        - ADMM dynamic tolerance update factor (default 1.e-12)
729 . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2)
730 . -tao_admm_minimum_spectral_penalty       - Set ADMM minimum spectral penalty (default 0)
731 . -tao_admm_dual_update                    - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic")
732 - -tao_admm_regularizer_type               - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold")
733 
734   Level: beginner
735 
736 .seealso: `TaoADMMSetMisfitHessianChangeStatus()`, `TaoADMMSetRegHessianChangeStatus()`, `TaoADMMGetSpectralPenalty()`,
737           `TaoADMMGetMisfitSubsolver()`, `TaoADMMGetRegularizationSubsolver()`, `TaoADMMSetConstraintVectorRHS()`,
738           `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetRegularizerCoefficient()`,
739           `TaoADMMSetRegularizerConstraintJacobian()`, `TaoADMMSetMisfitConstraintJacobian()`,
740           `TaoADMMSetMisfitObjectiveAndGradientRoutine()`, `TaoADMMSetMisfitHessianRoutine()`,
741           `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, `TaoADMMSetRegularizerHessianRoutine()`,
742           `TaoGetADMMParentTao()`, `TaoADMMGetDualVector()`, `TaoADMMSetRegularizerType()`,
743           `TaoADMMGetRegularizerType()`, `TaoADMMSetUpdateType()`, `TaoADMMGetUpdateType()`
744 M*/
745 
746 PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao)
747 {
748   TAO_ADMM       *am;
749 
750   PetscFunctionBegin;
751   PetscCall(PetscNewLog(tao,&am));
752 
753   tao->ops->destroy        = TaoDestroy_ADMM;
754   tao->ops->setup          = TaoSetUp_ADMM;
755   tao->ops->setfromoptions = TaoSetFromOptions_ADMM;
756   tao->ops->view           = TaoView_ADMM;
757   tao->ops->solve          = TaoSolve_ADMM;
758 
759   tao->data           = (void*)am;
760   am->l1epsilon       = 1e-6;
761   am->lambda          = 1e-4;
762   am->mu              = 1.;
763   am->muold           = 0.;
764   am->mueps           = PETSC_MACHINE_EPSILON;
765   am->mumin           = 0.;
766   am->orthval         = 0.2;
767   am->T               = 2;
768   am->parent          = tao;
769   am->update          = TAO_ADMM_UPDATE_BASIC;
770   am->regswitch       = TAO_ADMM_REGULARIZER_SOFT_THRESH;
771   am->tol             = PETSC_SMALL;
772   am->const_norm      = 0;
773   am->resnorm         = 0;
774   am->dualres         = 0;
775   am->ops->regobjgrad = NULL;
776   am->ops->reghess    = NULL;
777   am->gamma           = 1;
778   am->regobjgradP     = NULL;
779   am->reghessP        = NULL;
780   am->gatol_admm      = 1e-8;
781   am->catol_admm      = 0;
782   am->Hxchange        = PETSC_TRUE;
783   am->Hzchange        = PETSC_TRUE;
784   am->Hzbool          = PETSC_TRUE;
785   am->Hxbool          = PETSC_TRUE;
786 
787   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverX));
788   PetscCall(TaoSetOptionsPrefix(am->subsolverX,"misfit_"));
789   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX,(PetscObject)tao,1));
790   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverZ));
791   PetscCall(TaoSetOptionsPrefix(am->subsolverZ,"reg_"));
792   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ,(PetscObject)tao,1));
793 
794   PetscCall(TaoSetType(am->subsolverX,TAONLS));
795   PetscCall(TaoSetType(am->subsolverZ,TAONLS));
796   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", (PetscObject) tao));
797   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", (PetscObject) tao));
798   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",TaoADMMSetRegularizerType_ADMM));
799   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",TaoADMMGetRegularizerType_ADMM));
800   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",TaoADMMSetUpdateType_ADMM));
801   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",TaoADMMGetUpdateType_ADMM));
802   PetscFunctionReturn(0);
803 }
804 
805 /*@
806   TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines  whether Hessian matrix of misfit subsolver changes with respect to input vector.
807 
808   Collective on Tao
809 
810   Input Parameters:
811 +  tao - the Tao solver context.
812 -  b - the Hessian matrix change status boolean, PETSC_FALSE  when the Hessian matrix does not change, TRUE otherwise.
813 
814   Level: advanced
815 
816 .seealso: `TAOADMM`
817 
818 @*/
819 PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b)
820 {
821   TAO_ADMM *am = (TAO_ADMM*)tao->data;
822 
823   PetscFunctionBegin;
824   am->Hxchange = b;
825   PetscFunctionReturn(0);
826 }
827 
828 /*@
829   TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector.
830 
831   Collective on Tao
832 
833   Input Parameters:
834 +  tao - the Tao solver context
835 -  b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise.
836 
837   Level: advanced
838 
839 .seealso: `TAOADMM`
840 
841 @*/
842 PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b)
843 {
844   TAO_ADMM *am = (TAO_ADMM*)tao->data;
845 
846   PetscFunctionBegin;
847   am->Hzchange = b;
848   PetscFunctionReturn(0);
849 }
850 
851 /*@
852   TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value
853 
854   Collective on Tao
855 
856   Input Parameters:
857 +  tao - the Tao solver context
858 -  mu - spectral penalty
859 
860   Level: advanced
861 
862 .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TAOADMM`
863 @*/
864 PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu)
865 {
866   TAO_ADMM *am = (TAO_ADMM*)tao->data;
867 
868   PetscFunctionBegin;
869   am->mu = mu;
870   PetscFunctionReturn(0);
871 }
872 
873 /*@
874   TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value
875 
876   Collective on Tao
877 
878   Input Parameter:
879 .  tao - the Tao solver context
880 
881   Output Parameter:
882 .  mu - spectral penalty
883 
884   Level: advanced
885 
886 .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetSpectralPenalty()`, `TAOADMM`
887 @*/
888 PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu)
889 {
890   TAO_ADMM *am = (TAO_ADMM*)tao->data;
891 
892   PetscFunctionBegin;
893   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
894   PetscValidRealPointer(mu,2);
895   *mu = am->mu;
896   PetscFunctionReturn(0);
897 }
898 
899 /*@
900   TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM
901 
902   Collective on Tao
903 
904   Input Parameter:
905 .  tao - the Tao solver context
906 
907    Output Parameter:
908 .  misfit - the Tao subsolver context
909 
910   Level: advanced
911 
912 .seealso: `TAOADMM`
913 
914 @*/
915 PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit)
916 {
917   TAO_ADMM *am = (TAO_ADMM*)tao->data;
918 
919   PetscFunctionBegin;
920   *misfit = am->subsolverX;
921   PetscFunctionReturn(0);
922 }
923 
924 /*@
925   TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM
926 
927   Collective on Tao
928 
929   Input Parameter:
930 .  tao - the Tao solver context
931 
932   Output Parameter:
933 .  reg - the Tao subsolver context
934 
935   Level: advanced
936 
937 .seealso: `TAOADMM`
938 
939 @*/
940 PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg)
941 {
942   TAO_ADMM *am = (TAO_ADMM*)tao->data;
943 
944   PetscFunctionBegin;
945   *reg = am->subsolverZ;
946   PetscFunctionReturn(0);
947 }
948 
949 /*@
950   TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM
951 
952   Collective on Tao
953 
954   Input Parameters:
955 + tao - the Tao solver context
956 - c - RHS vector
957 
958   Level: advanced
959 
960 .seealso: `TAOADMM`
961 
962 @*/
963 PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c)
964 {
965   TAO_ADMM *am = (TAO_ADMM*)tao->data;
966 
967   PetscFunctionBegin;
968   am->constraint = c;
969   PetscFunctionReturn(0);
970 }
971 
972 /*@
973   TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty
974 
975   Collective on Tao
976 
977   Input Parameters:
978 +  tao - the Tao solver context
979 -  mu  - minimum spectral penalty value
980 
981   Level: advanced
982 
983 .seealso: `TaoADMMGetSpectralPenalty()`, `TAOADMM`
984 @*/
985 PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu)
986 {
987   TAO_ADMM *am = (TAO_ADMM*)tao->data;
988 
989   PetscFunctionBegin;
990   am->mumin= mu;
991   PetscFunctionReturn(0);
992 }
993 
994 /*@
995   TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case
996 
997   Collective on Tao
998 
999   Input Parameters:
1000 +  tao - the Tao solver context
1001 -  lambda - L1-norm regularizer coefficient
1002 
1003   Level: advanced
1004 
1005 .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
1006 
1007 @*/
1008 PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda)
1009 {
1010   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1011 
1012   PetscFunctionBegin;
1013   am->lambda = lambda;
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 /*@C
1018   TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable.
1019 
1020   Collective on Tao
1021 
1022   Input Parameters:
1023 + tao - the Tao solver context
1024 . J - user-created regularizer constraint Jacobian matrix
1025 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix
1026 . func - function pointer for the regularizer constraint Jacobian update function
1027 - ctx - user context for the regularizer Hessian
1028 
1029   Level: advanced
1030 
1031 .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
1032 
1033 @*/
1034 PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1035 {
1036   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1037 
1038   PetscFunctionBegin;
1039   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1040   if (J) {
1041     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
1042     PetscCheckSameComm(tao,1,J,2);
1043   }
1044   if (Jpre) {
1045     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
1046     PetscCheckSameComm(tao,1,Jpre,3);
1047   }
1048   if (ctx)  am->misfitjacobianP = ctx;
1049   if (func) am->ops->misfitjac  = func;
1050 
1051   if (J) {
1052     PetscCall(PetscObjectReference((PetscObject)J));
1053     PetscCall(MatDestroy(&am->JA));
1054     am->JA = J;
1055   }
1056   if (Jpre) {
1057     PetscCall(PetscObjectReference((PetscObject)Jpre));
1058     PetscCall(MatDestroy(&am->JApre));
1059     am->JApre = Jpre;
1060   }
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 /*@C
1065   TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable.
1066 
1067   Collective on Tao
1068 
1069   Input Parameters:
1070 + tao - the Tao solver context
1071 . J - user-created regularizer constraint Jacobian matrix
1072 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix
1073 . func - function pointer for the regularizer constraint Jacobian update function
1074 - ctx - user context for the regularizer Hessian
1075 
1076   Level: advanced
1077 
1078 .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetMisfitConstraintJacobian()`, `TAOADMM`
1079 
1080 @*/
1081 PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1082 {
1083   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1084 
1085   PetscFunctionBegin;
1086   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1087   if (J) {
1088     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
1089     PetscCheckSameComm(tao,1,J,2);
1090   }
1091   if (Jpre) {
1092     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
1093     PetscCheckSameComm(tao,1,Jpre,3);
1094   }
1095   if (ctx)  am->regjacobianP = ctx;
1096   if (func) am->ops->regjac  = func;
1097 
1098   if (J) {
1099     PetscCall(PetscObjectReference((PetscObject)J));
1100     PetscCall(MatDestroy(&am->JB));
1101     am->JB = J;
1102   }
1103   if (Jpre) {
1104     PetscCall(PetscObjectReference((PetscObject)Jpre));
1105     PetscCall(MatDestroy(&am->JBpre));
1106     am->JBpre = Jpre;
1107   }
1108   PetscFunctionReturn(0);
1109 }
1110 
1111 /*@C
1112    TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function
1113 
1114    Collective on tao
1115 
1116    Input Parameters:
1117 +    tao - the Tao context
1118 .    func - function pointer for the misfit value and gradient evaluation
1119 -    ctx - user context for the misfit
1120 
1121    Level: advanced
1122 
1123 .seealso: `TAOADMM`
1124 
1125 @*/
1126 PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx)
1127 {
1128   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1129 
1130   PetscFunctionBegin;
1131   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1132   am->misfitobjgradP     = ctx;
1133   am->ops->misfitobjgrad = func;
1134   PetscFunctionReturn(0);
1135 }
1136 
1137 /*@C
1138    TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back
1139    function into the algorithm, to be used for subsolverX.
1140 
1141    Collective on tao
1142 
1143    Input Parameters:
1144 + tao - the Tao context
1145 . H - user-created matrix for the Hessian of the misfit term
1146 . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term
1147 . func - function pointer for the misfit Hessian evaluation
1148 - ctx - user context for the misfit Hessian
1149 
1150    Level: advanced
1151 
1152 .seealso: `TAOADMM`
1153 
1154 @*/
1155 PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1156 {
1157   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1158 
1159   PetscFunctionBegin;
1160   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1161   if (H) {
1162     PetscValidHeaderSpecific(H,MAT_CLASSID,2);
1163     PetscCheckSameComm(tao,1,H,2);
1164   }
1165   if (Hpre) {
1166     PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
1167     PetscCheckSameComm(tao,1,Hpre,3);
1168   }
1169   if (ctx) {
1170     am->misfithessP = ctx;
1171   }
1172   if (func) {
1173     am->ops->misfithess = func;
1174   }
1175   if (H) {
1176     PetscCall(PetscObjectReference((PetscObject)H));
1177     PetscCall(MatDestroy(&am->Hx));
1178     am->Hx = H;
1179   }
1180   if (Hpre) {
1181     PetscCall(PetscObjectReference((PetscObject)Hpre));
1182     PetscCall(MatDestroy(&am->Hxpre));
1183     am->Hxpre = Hpre;
1184   }
1185   PetscFunctionReturn(0);
1186 }
1187 
1188 /*@C
1189    TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function
1190 
1191    Collective on tao
1192 
1193    Input Parameters:
1194 + tao - the Tao context
1195 . func - function pointer for the regularizer value and gradient evaluation
1196 - ctx - user context for the regularizer
1197 
1198    Level: advanced
1199 
1200 .seealso: `TAOADMM`
1201 
1202 @*/
1203 PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx)
1204 {
1205   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1206 
1207   PetscFunctionBegin;
1208   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1209   am->regobjgradP     = ctx;
1210   am->ops->regobjgrad = func;
1211   PetscFunctionReturn(0);
1212 }
1213 
1214 /*@C
1215    TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back
1216    function, to be used for subsolverZ.
1217 
1218    Collective on tao
1219 
1220    Input Parameters:
1221 + tao - the Tao context
1222 . H - user-created matrix for the Hessian of the regularization term
1223 . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term
1224 . func - function pointer for the regularizer Hessian evaluation
1225 - ctx - user context for the regularizer Hessian
1226 
1227    Level: advanced
1228 
1229 .seealso: `TAOADMM`
1230 
1231 @*/
1232 PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1233 {
1234   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1235 
1236   PetscFunctionBegin;
1237   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1238   if (H) {
1239     PetscValidHeaderSpecific(H,MAT_CLASSID,2);
1240     PetscCheckSameComm(tao,1,H,2);
1241   }
1242   if (Hpre) {
1243     PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
1244     PetscCheckSameComm(tao,1,Hpre,3);
1245   }
1246   if (ctx) {
1247     am->reghessP = ctx;
1248   }
1249   if (func) {
1250     am->ops->reghess = func;
1251   }
1252   if (H) {
1253     PetscCall(PetscObjectReference((PetscObject)H));
1254     PetscCall(MatDestroy(&am->Hz));
1255     am->Hz = H;
1256   }
1257   if (Hpre) {
1258     PetscCall(PetscObjectReference((PetscObject)Hpre));
1259     PetscCall(MatDestroy(&am->Hzpre));
1260     am->Hzpre = Hpre;
1261   }
1262   PetscFunctionReturn(0);
1263 }
1264 
1265 /*@
1266    TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver.
1267 
1268    Collective on tao
1269 
1270    Input Parameter:
1271 . tao - the Tao context
1272 
1273    Output Parameter:
1274 . admm_tao - the parent Tao context
1275 
1276    Level: advanced
1277 
1278 .seealso: `TAOADMM`
1279 
1280 @*/
1281 PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao)
1282 {
1283   PetscFunctionBegin;
1284   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1285   PetscCall(PetscObjectQuery((PetscObject)tao,"TaoGetADMMParentTao_ADMM", (PetscObject*) admm_tao));
1286   PetscFunctionReturn(0);
1287 }
1288 
1289 /*@
1290   TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state
1291 
1292   Not Collective
1293 
1294   Input Parameter:
1295 . tao - the Tao context
1296 
1297   Output Parameter:
1298 . Y - the current solution
1299 
1300   Level: intermediate
1301 
1302 .seealso: `TAOADMM`
1303 
1304 @*/
1305 PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y)
1306 {
1307   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1308 
1309   PetscFunctionBegin;
1310   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1311   *Y = am->y;
1312   PetscFunctionReturn(0);
1313 }
1314 
1315 /*@
1316   TaoADMMSetRegularizerType - Set regularizer type for ADMM routine
1317 
1318   Not Collective
1319 
1320   Input Parameters:
1321 + tao  - the Tao context
1322 - type - regularizer type
1323 
1324   Options Database:
1325 .  -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer
1326 
1327   Level: intermediate
1328 
1329 .seealso: `TaoADMMGetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
1330 @*/
1331 PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type)
1332 {
1333   PetscFunctionBegin;
1334   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1335   PetscValidLogicalCollectiveEnum(tao,type,2);
1336   PetscTryMethod(tao,"TaoADMMSetRegularizerType_C",(Tao,TaoADMMRegularizerType),(tao,type));
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 /*@
1341    TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM
1342 
1343    Not Collective
1344 
1345    Input Parameter:
1346 .  tao - the Tao context
1347 
1348    Output Parameter:
1349 .  type - the type of regularizer
1350 
1351    Level: intermediate
1352 
1353 .seealso: `TaoADMMSetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
1354 @*/
1355 PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type)
1356 {
1357   PetscFunctionBegin;
1358   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1359   PetscUseMethod(tao,"TaoADMMGetRegularizerType_C",(Tao,TaoADMMRegularizerType*),(tao,type));
1360   PetscFunctionReturn(0);
1361 }
1362 
1363 /*@
1364   TaoADMMSetUpdateType - Set update routine for ADMM routine
1365 
1366   Not Collective
1367 
1368   Input Parameters:
1369 + tao  - the Tao context
1370 - type - spectral parameter update type
1371 
1372   Level: intermediate
1373 
1374 .seealso: `TaoADMMGetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
1375 @*/
1376 PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type)
1377 {
1378   PetscFunctionBegin;
1379   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1380   PetscValidLogicalCollectiveEnum(tao,type,2);
1381   PetscTryMethod(tao,"TaoADMMSetUpdateType_C",(Tao,TaoADMMUpdateType),(tao,type));
1382   PetscFunctionReturn(0);
1383 }
1384 
1385 /*@
1386    TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM
1387 
1388    Not Collective
1389 
1390    Input Parameter:
1391 .  tao - the Tao context
1392 
1393    Output Parameter:
1394 .  type - the type of spectral penalty update routine
1395 
1396    Level: intermediate
1397 
1398 .seealso: `TaoADMMSetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
1399 @*/
1400 PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type)
1401 {
1402   PetscFunctionBegin;
1403   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1404   PetscUseMethod(tao,"TaoADMMGetUpdateType_C",(Tao,TaoADMMUpdateType*),(tao,type));
1405   PetscFunctionReturn(0);
1406 }
1407