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