xref: /petsc/src/tao/constrained/impls/admm/admm.c (revision 87edede185c0c71d403c393b2a06c8dc907b3f11)
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       break;
384     case (TAO_ADMM_REGULARIZER_SOFT_THRESH):
385       /* Soft Threshold. */
386       break;
387     }
388     if (am->ops->regobjgrad) PetscCall(TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao));
389     if (am->Hz && am->ops->reghess) {
390       PetscCall(TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao));
391     }
392   }
393 
394   switch (am->update) {
395   case TAO_ADMM_UPDATE_BASIC:
396     if (am->subsolverX->hessian) {
397       /* In basic case, Hessian does not get updated w.r.t. to spectral penalty
398        * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/
399       if (!am->xJI) {
400         PetscCall(MatAXPY(am->subsolverX->hessian,am->mu,am->ATA,DIFFERENT_NONZERO_PATTERN));
401       } else {
402         PetscCall(MatShift(am->subsolverX->hessian,am->mu));
403       }
404     }
405     if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) {
406       if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) {
407         PetscCall(MatAXPY(am->subsolverZ->hessian,am->mu,am->BTB,DIFFERENT_NONZERO_PATTERN));
408       } else {
409         PetscCall(MatShift(am->subsolverZ->hessian,am->mu));
410       }
411     }
412     break;
413   case TAO_ADMM_UPDATE_ADAPTIVE:
414   case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
415     break;
416   }
417 
418   PetscCall(PetscCitationsRegister(citation,&cited));
419   tao->reason = TAO_CONTINUE_ITERATING;
420 
421   while (tao->reason == TAO_CONTINUE_ITERATING) {
422     PetscTryTypeMethod(tao,update, tao->niter, tao->user_update);
423     PetscCall(VecCopy(am->Bz, am->Bzold));
424 
425     /* x update */
426     PetscCall(TaoSolve(am->subsolverX));
427     PetscCall(TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre));
428     PetscCall(MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution,am->Ax));
429 
430     am->Hxbool = PETSC_TRUE;
431 
432     /* z update */
433     switch (am->regswitch) {
434     case TAO_ADMM_REGULARIZER_USER:
435       PetscCall(TaoSolve(am->subsolverZ));
436       break;
437     case TAO_ADMM_REGULARIZER_SOFT_THRESH:
438       /* L1 assumes A,B jacobians are identity nxn matrix */
439       PetscCall(VecWAXPY(am->workJacobianRight,1/am->mu,am->y,am->Ax));
440       PetscCall(TaoSoftThreshold(am->workJacobianRight,-am->lambda/am->mu,am->lambda/am->mu,am->subsolverZ->solution));
441       break;
442     }
443     am->Hzbool = PETSC_TRUE;
444     /* Returns Ax + Bz - c with updated Ax,Bz vectors */
445     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
446     /* Dual variable, y += y + mu*(Ax+Bz-c) */
447     PetscCall(VecWAXPY(am->y, am->mu, am->residual, am->yold));
448 
449     /* stopping tolerance update */
450     PetscCall(TaoADMMToleranceUpdate(tao));
451 
452     /* Updating Spectral Penalty */
453     switch (am->update) {
454     case TAO_ADMM_UPDATE_BASIC:
455       am->muold = am->mu;
456       break;
457     case TAO_ADMM_UPDATE_ADAPTIVE:
458     case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED:
459       if (tao->niter == 0) {
460         PetscCall(VecCopy(am->y, am->y0));
461         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
462         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
463         PetscCall(VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold));
464         PetscCall(VecCopy(am->Ax, am->Axold));
465         PetscCall(VecCopy(am->Bz, am->Bz0));
466         am->muold = am->mu;
467       } else if (tao->niter % am->T == 1) {
468         /* we have compute Bzold in a previous iteration, and we computed Ax above */
469         PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold));
470         if (am->constraint) PetscCall(VecAXPY(am->residual, -1., am->constraint));
471         PetscCall(VecWAXPY(am->yhat, -am->mu, am->residual, am->yold));
472         PetscCall(AdaptiveADMMPenaltyUpdate(tao));
473         PetscCall(VecCopy(am->Ax, am->Axold));
474         PetscCall(VecCopy(am->Bz, am->Bz0));
475         PetscCall(VecCopy(am->yhat, am->yhatold));
476         PetscCall(VecCopy(am->y, am->y0));
477       } else {
478         am->muold = am->mu;
479       }
480       break;
481     default:
482       break;
483     }
484     tao->niter++;
485 
486     /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/
487     switch (am->regswitch) {
488     case TAO_ADMM_REGULARIZER_USER:
489       if (is_reg_shell) {
490         PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,&reg_func));
491       } else {
492         (*am->ops->regobjgrad)(am->subsolverZ,am->subsolverX->solution,&reg_func,tempL,am->regobjgradP);
493       }
494       break;
495     case TAO_ADMM_REGULARIZER_SOFT_THRESH:
496       PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,&reg_func));
497       break;
498     }
499     PetscCall(VecCopy(am->y,am->yold));
500     PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual));
501     PetscCall(VecNorm(am->residual,NORM_2,&am->resnorm));
502     PetscCall(TaoLogConvergenceHistory(tao,am->last_misfit_val + reg_func,am->dualres,am->resnorm,tao->ksp_its));
503 
504     PetscCall(TaoMonitor(tao,tao->niter,am->last_misfit_val + reg_func,am->dualres,am->resnorm,1.0));
505     PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
506   }
507   /* Update vectors */
508   PetscCall(VecCopy(am->subsolverX->solution,tao->solution));
509   PetscCall(VecCopy(am->subsolverX->gradient,tao->gradient));
510   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", NULL));
511   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", NULL));
512   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL));
513   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL));
514   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL));
515   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL));
516   PetscFunctionReturn(0);
517 }
518 
519 static PetscErrorCode TaoSetFromOptions_ADMM(Tao tao,PetscOptionItems *PetscOptionsObject)
520 {
521   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
522 
523   PetscFunctionBegin;
524   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. ");
525   PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient","regularizer constant","",am->lambda,&am->lambda,NULL));
526   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty","Constant for Augmented Lagrangian term.","",am->mu,&am->mu,NULL));
527   PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter","x relaxation parameter for Z update.","",am->gamma,&am->gamma,NULL));
528   PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor","ADMM dynamic tolerance update factor.","",am->tol,&am->tol,NULL));
529   PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor","ADMM spectral penalty update curvature safeguard value.","",am->orthval,&am->orthval,NULL));
530   PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty","Set ADMM minimum spectral penalty.","",am->mumin,&am->mumin,NULL));
531   PetscCall(PetscOptionsEnum("-tao_admm_dual_update","Lagrangian dual update policy","TaoADMMUpdateType",TaoADMMUpdateTypes,(PetscEnum)am->update,(PetscEnum*)&am->update,NULL));
532   PetscCall(PetscOptionsEnum("-tao_admm_regularizer_type","ADMM regularizer update rule","TaoADMMRegularizerType",TaoADMMRegularizerTypes,(PetscEnum)am->regswitch,(PetscEnum*)&am->regswitch,NULL));
533   PetscOptionsHeadEnd();
534   PetscCall(TaoSetFromOptions(am->subsolverX));
535   if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) {
536     PetscCall(TaoSetFromOptions(am->subsolverZ));
537   }
538   PetscFunctionReturn(0);
539 }
540 
541 static PetscErrorCode TaoView_ADMM(Tao tao,PetscViewer viewer)
542 {
543   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
544 
545   PetscFunctionBegin;
546   PetscCall(PetscViewerASCIIPushTab(viewer));
547   PetscCall(TaoView(am->subsolverX,viewer));
548   PetscCall(TaoView(am->subsolverZ,viewer));
549   PetscCall(PetscViewerASCIIPopTab(viewer));
550   PetscFunctionReturn(0);
551 }
552 
553 static PetscErrorCode TaoSetUp_ADMM(Tao tao)
554 {
555   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
556   PetscInt       n,N,M;
557 
558   PetscFunctionBegin;
559   PetscCall(VecGetLocalSize(tao->solution,&n));
560   PetscCall(VecGetSize(tao->solution,&N));
561   /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */
562   if (!am->JB) {
563     am->zJI   = PETSC_TRUE;
564     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JB));
565     PetscCall(MatShellSetOperation(am->JB,MATOP_MULT,(void (*)(void))JacobianIdentityB));
566     PetscCall(MatShellSetOperation(am->JB,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentityB));
567     am->JBpre = am->JB;
568   }
569   if (!am->JA) {
570     am->xJI   = PETSC_TRUE;
571     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JA));
572     PetscCall(MatShellSetOperation(am->JA,MATOP_MULT,(void (*)(void))JacobianIdentity));
573     PetscCall(MatShellSetOperation(am->JA,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentity));
574     am->JApre = am->JA;
575   }
576   PetscCall(MatCreateVecs(am->JA,NULL,&am->Ax));
577   if (!tao->gradient) {
578     PetscCall(VecDuplicate(tao->solution,&tao->gradient));
579   }
580   PetscCall(TaoSetSolution(am->subsolverX, tao->solution));
581   if (!am->z) {
582     PetscCall(VecDuplicate(tao->solution,&am->z));
583     PetscCall(VecSet(am->z,0.0));
584   }
585   PetscCall(TaoSetSolution(am->subsolverZ, am->z));
586   if (!am->workLeft) {
587     PetscCall(VecDuplicate(tao->solution,&am->workLeft));
588   }
589   if (!am->Axold) {
590     PetscCall(VecDuplicate(am->Ax,&am->Axold));
591   }
592   if (!am->workJacobianRight) {
593     PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight));
594   }
595   if (!am->workJacobianRight2) {
596     PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight2));
597   }
598   if (!am->Bz) {
599     PetscCall(VecDuplicate(am->Ax,&am->Bz));
600   }
601   if (!am->Bzold) {
602     PetscCall(VecDuplicate(am->Ax,&am->Bzold));
603   }
604   if (!am->Bz0) {
605     PetscCall(VecDuplicate(am->Ax,&am->Bz0));
606   }
607   if (!am->y) {
608     PetscCall(VecDuplicate(am->Ax,&am->y));
609     PetscCall(VecSet(am->y,0.0));
610   }
611   if (!am->yold) {
612     PetscCall(VecDuplicate(am->Ax,&am->yold));
613     PetscCall(VecSet(am->yold,0.0));
614   }
615   if (!am->y0) {
616     PetscCall(VecDuplicate(am->Ax,&am->y0));
617     PetscCall(VecSet(am->y0,0.0));
618   }
619   if (!am->yhat) {
620     PetscCall(VecDuplicate(am->Ax,&am->yhat));
621     PetscCall(VecSet(am->yhat,0.0));
622   }
623   if (!am->yhatold) {
624     PetscCall(VecDuplicate(am->Ax,&am->yhatold));
625     PetscCall(VecSet(am->yhatold,0.0));
626   }
627   if (!am->residual) {
628     PetscCall(VecDuplicate(am->Ax,&am->residual));
629     PetscCall(VecSet(am->residual,0.0));
630   }
631   if (!am->constraint) {
632     am->constraint = NULL;
633   } else {
634     PetscCall(VecGetSize(am->constraint,&M));
635     PetscCheck(M == N,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Solution vector and constraint vector must be of same size!");
636   }
637 
638   /* Save changed tao tolerance for adaptive tolerance */
639   if (tao->gatol_changed) {
640     am->gatol_admm = tao->gatol;
641   }
642   if (tao->catol_changed) {
643     am->catol_admm = tao->catol;
644   }
645 
646   /*Update spectral and dual elements to X subsolver */
647   PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao));
648   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX,am->JA,am->JApre, am->ops->misfitjac, am->misfitjacobianP));
649   PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ,am->JB,am->JBpre, am->ops->regjac, am->regjacobianP));
650   PetscFunctionReturn(0);
651 }
652 
653 static PetscErrorCode TaoDestroy_ADMM(Tao tao)
654 {
655   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
656 
657   PetscFunctionBegin;
658   PetscCall(VecDestroy(&am->z));
659   PetscCall(VecDestroy(&am->Ax));
660   PetscCall(VecDestroy(&am->Axold));
661   PetscCall(VecDestroy(&am->Bz));
662   PetscCall(VecDestroy(&am->Bzold));
663   PetscCall(VecDestroy(&am->Bz0));
664   PetscCall(VecDestroy(&am->residual));
665   PetscCall(VecDestroy(&am->y));
666   PetscCall(VecDestroy(&am->yold));
667   PetscCall(VecDestroy(&am->y0));
668   PetscCall(VecDestroy(&am->yhat));
669   PetscCall(VecDestroy(&am->yhatold));
670   PetscCall(VecDestroy(&am->workLeft));
671   PetscCall(VecDestroy(&am->workJacobianRight));
672   PetscCall(VecDestroy(&am->workJacobianRight2));
673 
674   PetscCall(MatDestroy(&am->JA));
675   PetscCall(MatDestroy(&am->JB));
676   if (!am->xJI) {
677     PetscCall(MatDestroy(&am->JApre));
678   }
679   if (!am->zJI) {
680     PetscCall(MatDestroy(&am->JBpre));
681   }
682   if (am->Hx) {
683     PetscCall(MatDestroy(&am->Hx));
684     PetscCall(MatDestroy(&am->Hxpre));
685   }
686   if (am->Hz) {
687     PetscCall(MatDestroy(&am->Hz));
688     PetscCall(MatDestroy(&am->Hzpre));
689   }
690   PetscCall(MatDestroy(&am->ATA));
691   PetscCall(MatDestroy(&am->BTB));
692   PetscCall(TaoDestroy(&am->subsolverX));
693   PetscCall(TaoDestroy(&am->subsolverZ));
694   am->parent = NULL;
695   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL));
696   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL));
697   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL));
698   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL));
699   PetscCall(PetscFree(tao->data));
700   PetscFunctionReturn(0);
701 }
702 
703 /*MC
704 
705   TAOADMM - Alternating direction method of multipliers method fo solving linear problems with
706             constraints. in a min_x f(x) + g(z)  s.t. Ax+Bz=c.
707             This algorithm employs two sub Tao solvers, of which type can be specified
708             by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers.
709             Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via
710             TaoADMMSet{Misfit,Regularizer}HessianChangeStatus.
711             Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type.
712             There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update,
713             currently there are basic option and adaptive option.
714             Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian.
715             c can be set with TaoADMMSetConstraintVectorRHS.
716             The user can also provide regularizer weight for second subsolver.
717 
718   References:
719 . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom
720           "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation"
721           The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017.
722 
723   Options Database Keys:
724 + -tao_admm_regularizer_coefficient        - regularizer constant (default 1.e-6)
725 . -tao_admm_spectral_penalty               - Constant for Augmented Lagrangian term (default 1.)
726 . -tao_admm_relaxation_parameter           - relaxation parameter for Z update (default 1.)
727 . -tao_admm_tolerance_update_factor        - ADMM dynamic tolerance update factor (default 1.e-12)
728 . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2)
729 . -tao_admm_minimum_spectral_penalty       - Set ADMM minimum spectral penalty (default 0)
730 . -tao_admm_dual_update                    - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic")
731 - -tao_admm_regularizer_type               - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold")
732 
733   Level: beginner
734 
735 .seealso: `TaoADMMSetMisfitHessianChangeStatus()`, `TaoADMMSetRegHessianChangeStatus()`, `TaoADMMGetSpectralPenalty()`,
736           `TaoADMMGetMisfitSubsolver()`, `TaoADMMGetRegularizationSubsolver()`, `TaoADMMSetConstraintVectorRHS()`,
737           `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetRegularizerCoefficient()`,
738           `TaoADMMSetRegularizerConstraintJacobian()`, `TaoADMMSetMisfitConstraintJacobian()`,
739           `TaoADMMSetMisfitObjectiveAndGradientRoutine()`, `TaoADMMSetMisfitHessianRoutine()`,
740           `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, `TaoADMMSetRegularizerHessianRoutine()`,
741           `TaoGetADMMParentTao()`, `TaoADMMGetDualVector()`, `TaoADMMSetRegularizerType()`,
742           `TaoADMMGetRegularizerType()`, `TaoADMMSetUpdateType()`, `TaoADMMGetUpdateType()`
743 M*/
744 
745 PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao)
746 {
747   TAO_ADMM       *am;
748 
749   PetscFunctionBegin;
750   PetscCall(PetscNewLog(tao,&am));
751 
752   tao->ops->destroy        = TaoDestroy_ADMM;
753   tao->ops->setup          = TaoSetUp_ADMM;
754   tao->ops->setfromoptions = TaoSetFromOptions_ADMM;
755   tao->ops->view           = TaoView_ADMM;
756   tao->ops->solve          = TaoSolve_ADMM;
757 
758   tao->data           = (void*)am;
759   am->l1epsilon       = 1e-6;
760   am->lambda          = 1e-4;
761   am->mu              = 1.;
762   am->muold           = 0.;
763   am->mueps           = PETSC_MACHINE_EPSILON;
764   am->mumin           = 0.;
765   am->orthval         = 0.2;
766   am->T               = 2;
767   am->parent          = tao;
768   am->update          = TAO_ADMM_UPDATE_BASIC;
769   am->regswitch       = TAO_ADMM_REGULARIZER_SOFT_THRESH;
770   am->tol             = PETSC_SMALL;
771   am->const_norm      = 0;
772   am->resnorm         = 0;
773   am->dualres         = 0;
774   am->ops->regobjgrad = NULL;
775   am->ops->reghess    = NULL;
776   am->gamma           = 1;
777   am->regobjgradP     = NULL;
778   am->reghessP        = NULL;
779   am->gatol_admm      = 1e-8;
780   am->catol_admm      = 0;
781   am->Hxchange        = PETSC_TRUE;
782   am->Hzchange        = PETSC_TRUE;
783   am->Hzbool          = PETSC_TRUE;
784   am->Hxbool          = PETSC_TRUE;
785 
786   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverX));
787   PetscCall(TaoSetOptionsPrefix(am->subsolverX,"misfit_"));
788   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX,(PetscObject)tao,1));
789   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverZ));
790   PetscCall(TaoSetOptionsPrefix(am->subsolverZ,"reg_"));
791   PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ,(PetscObject)tao,1));
792 
793   PetscCall(TaoSetType(am->subsolverX,TAONLS));
794   PetscCall(TaoSetType(am->subsolverZ,TAONLS));
795   PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", (PetscObject) tao));
796   PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", (PetscObject) tao));
797   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",TaoADMMSetRegularizerType_ADMM));
798   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",TaoADMMGetRegularizerType_ADMM));
799   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",TaoADMMSetUpdateType_ADMM));
800   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",TaoADMMGetUpdateType_ADMM));
801   PetscFunctionReturn(0);
802 }
803 
804 /*@
805   TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines  whether Hessian matrix of misfit subsolver changes with respect to input vector.
806 
807   Collective on Tao
808 
809   Input Parameters:
810 +  tao - the Tao solver context.
811 -  b - the Hessian matrix change status boolean, PETSC_FALSE  when the Hessian matrix does not change, TRUE otherwise.
812 
813   Level: advanced
814 
815 .seealso: `TAOADMM`
816 
817 @*/
818 PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b)
819 {
820   TAO_ADMM *am = (TAO_ADMM*)tao->data;
821 
822   PetscFunctionBegin;
823   am->Hxchange = b;
824   PetscFunctionReturn(0);
825 }
826 
827 /*@
828   TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector.
829 
830   Collective on Tao
831 
832   Input Parameters:
833 +  tao - the Tao solver context
834 -  b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise.
835 
836   Level: advanced
837 
838 .seealso: `TAOADMM`
839 
840 @*/
841 PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b)
842 {
843   TAO_ADMM *am = (TAO_ADMM*)tao->data;
844 
845   PetscFunctionBegin;
846   am->Hzchange = b;
847   PetscFunctionReturn(0);
848 }
849 
850 /*@
851   TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value
852 
853   Collective on Tao
854 
855   Input Parameters:
856 +  tao - the Tao solver context
857 -  mu - spectral penalty
858 
859   Level: advanced
860 
861 .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TAOADMM`
862 @*/
863 PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu)
864 {
865   TAO_ADMM *am = (TAO_ADMM*)tao->data;
866 
867   PetscFunctionBegin;
868   am->mu = mu;
869   PetscFunctionReturn(0);
870 }
871 
872 /*@
873   TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value
874 
875   Collective on Tao
876 
877   Input Parameter:
878 .  tao - the Tao solver context
879 
880   Output Parameter:
881 .  mu - spectral penalty
882 
883   Level: advanced
884 
885 .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetSpectralPenalty()`, `TAOADMM`
886 @*/
887 PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu)
888 {
889   TAO_ADMM *am = (TAO_ADMM*)tao->data;
890 
891   PetscFunctionBegin;
892   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
893   PetscValidRealPointer(mu,2);
894   *mu = am->mu;
895   PetscFunctionReturn(0);
896 }
897 
898 /*@
899   TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM
900 
901   Collective on Tao
902 
903   Input Parameter:
904 .  tao - the Tao solver context
905 
906    Output Parameter:
907 .  misfit - the Tao subsolver context
908 
909   Level: advanced
910 
911 .seealso: `TAOADMM`
912 
913 @*/
914 PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit)
915 {
916   TAO_ADMM *am = (TAO_ADMM*)tao->data;
917 
918   PetscFunctionBegin;
919   *misfit = am->subsolverX;
920   PetscFunctionReturn(0);
921 }
922 
923 /*@
924   TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM
925 
926   Collective on Tao
927 
928   Input Parameter:
929 .  tao - the Tao solver context
930 
931   Output Parameter:
932 .  reg - the Tao subsolver context
933 
934   Level: advanced
935 
936 .seealso: `TAOADMM`
937 
938 @*/
939 PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg)
940 {
941   TAO_ADMM *am = (TAO_ADMM*)tao->data;
942 
943   PetscFunctionBegin;
944   *reg = am->subsolverZ;
945   PetscFunctionReturn(0);
946 }
947 
948 /*@
949   TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM
950 
951   Collective on Tao
952 
953   Input Parameters:
954 + tao - the Tao solver context
955 - c - RHS vector
956 
957   Level: advanced
958 
959 .seealso: `TAOADMM`
960 
961 @*/
962 PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c)
963 {
964   TAO_ADMM *am = (TAO_ADMM*)tao->data;
965 
966   PetscFunctionBegin;
967   am->constraint = c;
968   PetscFunctionReturn(0);
969 }
970 
971 /*@
972   TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty
973 
974   Collective on Tao
975 
976   Input Parameters:
977 +  tao - the Tao solver context
978 -  mu  - minimum spectral penalty value
979 
980   Level: advanced
981 
982 .seealso: `TaoADMMGetSpectralPenalty()`, `TAOADMM`
983 @*/
984 PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu)
985 {
986   TAO_ADMM *am = (TAO_ADMM*)tao->data;
987 
988   PetscFunctionBegin;
989   am->mumin= mu;
990   PetscFunctionReturn(0);
991 }
992 
993 /*@
994   TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case
995 
996   Collective on Tao
997 
998   Input Parameters:
999 +  tao - the Tao solver context
1000 -  lambda - L1-norm regularizer coefficient
1001 
1002   Level: advanced
1003 
1004 .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
1005 
1006 @*/
1007 PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda)
1008 {
1009   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1010 
1011   PetscFunctionBegin;
1012   am->lambda = lambda;
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 /*@C
1017   TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable.
1018 
1019   Collective on Tao
1020 
1021   Input Parameters:
1022 + tao - the Tao solver context
1023 . J - user-created regularizer constraint Jacobian matrix
1024 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix
1025 . func - function pointer for the regularizer constraint Jacobian update function
1026 - ctx - user context for the regularizer Hessian
1027 
1028   Level: advanced
1029 
1030 .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM`
1031 
1032 @*/
1033 PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1034 {
1035   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1036 
1037   PetscFunctionBegin;
1038   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1039   if (J) {
1040     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
1041     PetscCheckSameComm(tao,1,J,2);
1042   }
1043   if (Jpre) {
1044     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
1045     PetscCheckSameComm(tao,1,Jpre,3);
1046   }
1047   if (ctx)  am->misfitjacobianP = ctx;
1048   if (func) am->ops->misfitjac  = func;
1049 
1050   if (J) {
1051     PetscCall(PetscObjectReference((PetscObject)J));
1052     PetscCall(MatDestroy(&am->JA));
1053     am->JA = J;
1054   }
1055   if (Jpre) {
1056     PetscCall(PetscObjectReference((PetscObject)Jpre));
1057     PetscCall(MatDestroy(&am->JApre));
1058     am->JApre = Jpre;
1059   }
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 /*@C
1064   TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable.
1065 
1066   Collective on Tao
1067 
1068   Input Parameters:
1069 + tao - the Tao solver context
1070 . J - user-created regularizer constraint Jacobian matrix
1071 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix
1072 . func - function pointer for the regularizer constraint Jacobian update function
1073 - ctx - user context for the regularizer Hessian
1074 
1075   Level: advanced
1076 
1077 .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetMisfitConstraintJacobian()`, `TAOADMM`
1078 
1079 @*/
1080 PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1081 {
1082   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1083 
1084   PetscFunctionBegin;
1085   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1086   if (J) {
1087     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
1088     PetscCheckSameComm(tao,1,J,2);
1089   }
1090   if (Jpre) {
1091     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
1092     PetscCheckSameComm(tao,1,Jpre,3);
1093   }
1094   if (ctx)  am->regjacobianP = ctx;
1095   if (func) am->ops->regjac  = func;
1096 
1097   if (J) {
1098     PetscCall(PetscObjectReference((PetscObject)J));
1099     PetscCall(MatDestroy(&am->JB));
1100     am->JB = J;
1101   }
1102   if (Jpre) {
1103     PetscCall(PetscObjectReference((PetscObject)Jpre));
1104     PetscCall(MatDestroy(&am->JBpre));
1105     am->JBpre = Jpre;
1106   }
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 /*@C
1111    TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function
1112 
1113    Collective on tao
1114 
1115    Input Parameters:
1116 +    tao - the Tao context
1117 .    func - function pointer for the misfit value and gradient evaluation
1118 -    ctx - user context for the misfit
1119 
1120    Level: advanced
1121 
1122 .seealso: `TAOADMM`
1123 
1124 @*/
1125 PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx)
1126 {
1127   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1128 
1129   PetscFunctionBegin;
1130   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1131   am->misfitobjgradP     = ctx;
1132   am->ops->misfitobjgrad = func;
1133   PetscFunctionReturn(0);
1134 }
1135 
1136 /*@C
1137    TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back
1138    function into the algorithm, to be used for subsolverX.
1139 
1140    Collective on tao
1141 
1142    Input Parameters:
1143 + tao - the Tao context
1144 . H - user-created matrix for the Hessian of the misfit term
1145 . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term
1146 . func - function pointer for the misfit Hessian evaluation
1147 - ctx - user context for the misfit Hessian
1148 
1149    Level: advanced
1150 
1151 .seealso: `TAOADMM`
1152 
1153 @*/
1154 PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1155 {
1156   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1157 
1158   PetscFunctionBegin;
1159   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1160   if (H) {
1161     PetscValidHeaderSpecific(H,MAT_CLASSID,2);
1162     PetscCheckSameComm(tao,1,H,2);
1163   }
1164   if (Hpre) {
1165     PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
1166     PetscCheckSameComm(tao,1,Hpre,3);
1167   }
1168   if (ctx) {
1169     am->misfithessP = ctx;
1170   }
1171   if (func) {
1172     am->ops->misfithess = func;
1173   }
1174   if (H) {
1175     PetscCall(PetscObjectReference((PetscObject)H));
1176     PetscCall(MatDestroy(&am->Hx));
1177     am->Hx = H;
1178   }
1179   if (Hpre) {
1180     PetscCall(PetscObjectReference((PetscObject)Hpre));
1181     PetscCall(MatDestroy(&am->Hxpre));
1182     am->Hxpre = Hpre;
1183   }
1184   PetscFunctionReturn(0);
1185 }
1186 
1187 /*@C
1188    TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function
1189 
1190    Collective on tao
1191 
1192    Input Parameters:
1193 + tao - the Tao context
1194 . func - function pointer for the regularizer value and gradient evaluation
1195 - ctx - user context for the regularizer
1196 
1197    Level: advanced
1198 
1199 .seealso: `TAOADMM`
1200 
1201 @*/
1202 PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx)
1203 {
1204   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1205 
1206   PetscFunctionBegin;
1207   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1208   am->regobjgradP     = ctx;
1209   am->ops->regobjgrad = func;
1210   PetscFunctionReturn(0);
1211 }
1212 
1213 /*@C
1214    TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back
1215    function, to be used for subsolverZ.
1216 
1217    Collective on tao
1218 
1219    Input Parameters:
1220 + tao - the Tao context
1221 . H - user-created matrix for the Hessian of the regularization term
1222 . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term
1223 . func - function pointer for the regularizer Hessian evaluation
1224 - ctx - user context for the regularizer Hessian
1225 
1226    Level: advanced
1227 
1228 .seealso: `TAOADMM`
1229 
1230 @*/
1231 PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
1232 {
1233   TAO_ADMM       *am = (TAO_ADMM*)tao->data;
1234 
1235   PetscFunctionBegin;
1236   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1237   if (H) {
1238     PetscValidHeaderSpecific(H,MAT_CLASSID,2);
1239     PetscCheckSameComm(tao,1,H,2);
1240   }
1241   if (Hpre) {
1242     PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
1243     PetscCheckSameComm(tao,1,Hpre,3);
1244   }
1245   if (ctx) {
1246     am->reghessP = ctx;
1247   }
1248   if (func) {
1249     am->ops->reghess = func;
1250   }
1251   if (H) {
1252     PetscCall(PetscObjectReference((PetscObject)H));
1253     PetscCall(MatDestroy(&am->Hz));
1254     am->Hz = H;
1255   }
1256   if (Hpre) {
1257     PetscCall(PetscObjectReference((PetscObject)Hpre));
1258     PetscCall(MatDestroy(&am->Hzpre));
1259     am->Hzpre = Hpre;
1260   }
1261   PetscFunctionReturn(0);
1262 }
1263 
1264 /*@
1265    TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver.
1266 
1267    Collective on tao
1268 
1269    Input Parameter:
1270 . tao - the Tao context
1271 
1272    Output Parameter:
1273 . admm_tao - the parent Tao context
1274 
1275    Level: advanced
1276 
1277 .seealso: `TAOADMM`
1278 
1279 @*/
1280 PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao)
1281 {
1282   PetscFunctionBegin;
1283   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1284   PetscCall(PetscObjectQuery((PetscObject)tao,"TaoGetADMMParentTao_ADMM", (PetscObject*) admm_tao));
1285   PetscFunctionReturn(0);
1286 }
1287 
1288 /*@
1289   TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state
1290 
1291   Not Collective
1292 
1293   Input Parameter:
1294 . tao - the Tao context
1295 
1296   Output Parameter:
1297 . Y - the current solution
1298 
1299   Level: intermediate
1300 
1301 .seealso: `TAOADMM`
1302 
1303 @*/
1304 PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y)
1305 {
1306   TAO_ADMM *am = (TAO_ADMM*)tao->data;
1307 
1308   PetscFunctionBegin;
1309   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1310   *Y = am->y;
1311   PetscFunctionReturn(0);
1312 }
1313 
1314 /*@
1315   TaoADMMSetRegularizerType - Set regularizer type for ADMM routine
1316 
1317   Not Collective
1318 
1319   Input Parameters:
1320 + tao  - the Tao context
1321 - type - regularizer type
1322 
1323   Options Database:
1324 .  -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer
1325 
1326   Level: intermediate
1327 
1328 .seealso: `TaoADMMGetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
1329 @*/
1330 PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type)
1331 {
1332   PetscFunctionBegin;
1333   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1334   PetscValidLogicalCollectiveEnum(tao,type,2);
1335   PetscTryMethod(tao,"TaoADMMSetRegularizerType_C",(Tao,TaoADMMRegularizerType),(tao,type));
1336   PetscFunctionReturn(0);
1337 }
1338 
1339 /*@
1340    TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM
1341 
1342    Not Collective
1343 
1344    Input Parameter:
1345 .  tao - the Tao context
1346 
1347    Output Parameter:
1348 .  type - the type of regularizer
1349 
1350    Level: intermediate
1351 
1352 .seealso: `TaoADMMSetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM`
1353 @*/
1354 PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type)
1355 {
1356   PetscFunctionBegin;
1357   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1358   PetscUseMethod(tao,"TaoADMMGetRegularizerType_C",(Tao,TaoADMMRegularizerType*),(tao,type));
1359   PetscFunctionReturn(0);
1360 }
1361 
1362 /*@
1363   TaoADMMSetUpdateType - Set update routine for ADMM routine
1364 
1365   Not Collective
1366 
1367   Input Parameters:
1368 + tao  - the Tao context
1369 - type - spectral parameter update type
1370 
1371   Level: intermediate
1372 
1373 .seealso: `TaoADMMGetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
1374 @*/
1375 PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type)
1376 {
1377   PetscFunctionBegin;
1378   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1379   PetscValidLogicalCollectiveEnum(tao,type,2);
1380   PetscTryMethod(tao,"TaoADMMSetUpdateType_C",(Tao,TaoADMMUpdateType),(tao,type));
1381   PetscFunctionReturn(0);
1382 }
1383 
1384 /*@
1385    TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM
1386 
1387    Not Collective
1388 
1389    Input Parameter:
1390 .  tao - the Tao context
1391 
1392    Output Parameter:
1393 .  type - the type of spectral penalty update routine
1394 
1395    Level: intermediate
1396 
1397 .seealso: `TaoADMMSetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM`
1398 @*/
1399 PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type)
1400 {
1401   PetscFunctionBegin;
1402   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1403   PetscUseMethod(tao,"TaoADMMGetUpdateType_C",(Tao,TaoADMMUpdateType*),(tao,type));
1404   PetscFunctionReturn(0);
1405 }
1406