Actual source code: admm.h
1: #pragma once
2: #include <petsc/private/taoimpl.h>
4: typedef struct _TaoADMMOps *TaoADMMOps;
6: struct _TaoADMMOps {
7: PetscErrorCode (*misfitobjgrad)(Tao, Vec, PetscReal *, Vec, void *);
8: PetscErrorCode (*misfithess)(Tao, Vec, Mat, Mat, void *);
9: PetscErrorCode (*misfitjac)(Tao, Vec, Mat, Mat, void *);
10: PetscErrorCode (*regobjgrad)(Tao, Vec, PetscReal *, Vec, void *);
11: PetscErrorCode (*reghess)(Tao, Vec, Mat, Mat, void *);
12: PetscErrorCode (*regjac)(Tao, Vec, Mat, Mat, void *);
13: };
15: typedef struct {
16: PETSCHEADER(struct _TaoADMMOps);
17: Tao subsolverX, subsolverZ, parent;
18: Vec residual, y, yold, y0, yhat, yhatold, constraint;
19: Vec z, zold, Ax, Bz, Axold, Bzold, Bz0;
20: Vec workLeft, workJacobianRight, workJacobianRight2; /*Ax,Bz,y,constraint are workJacobianRight sized. workLeft is solution sized */
21: Mat Hx, Hxpre, Hz, Hzpre, ATA, BTB, JA, JApre, JB, JBpre;
22: void *regobjgradP;
23: void *reghessP;
24: void *regjacobianP;
25: void *misfitobjgradP;
26: void *misfithessP;
27: void *misfitjacobianP;
28: PetscReal gamma, last_misfit_val, last_reg_val, l1epsilon;
29: PetscReal lambda, mu, muold, orthval, mueps, tol, mumin;
30: PetscReal Bzdiffnorm, dualres, resnorm, const_norm;
31: PetscReal gatol_admm, catol_admm;
32: PetscInt T; /* adaptive iteration cutoff */
33: PetscBool xJI, zJI; /* Bool to check whether A,B Jacobians are NULL-set identity */
34: PetscBool Hxchange, Hzchange; /* Bool to check whether Hx,Hz change wrt to x and z */
35: PetscBool Hxbool, Hzbool; /* Bool to make sure Hessian gets updated only once for Hchange False case */
36: TaoADMMUpdateType update; /* update policy for mu */
37: TaoADMMRegularizerType regswitch; /* regularization policy */
38: } TAO_ADMM;