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