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