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