xref: /petsc/src/tao/constrained/impls/admm/admm.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
1*a4963045SJacob Faibussowitsch #pragma once
26285c0a3SHansol  Suh #include <petsc/private/taoimpl.h>
36285c0a3SHansol  Suh 
46285c0a3SHansol  Suh typedef struct _TaoADMMOps *TaoADMMOps;
56285c0a3SHansol  Suh 
66285c0a3SHansol  Suh struct _TaoADMMOps {
76285c0a3SHansol  Suh   PetscErrorCode (*misfitobjgrad)(Tao, Vec, PetscReal *, Vec, void *);
86285c0a3SHansol  Suh   PetscErrorCode (*misfithess)(Tao, Vec, Mat, Mat, void *);
96285c0a3SHansol  Suh   PetscErrorCode (*misfitjac)(Tao, Vec, Mat, Mat, void *);
106285c0a3SHansol  Suh   PetscErrorCode (*regobjgrad)(Tao, Vec, PetscReal *, Vec, void *);
116285c0a3SHansol  Suh   PetscErrorCode (*reghess)(Tao, Vec, Mat, Mat, void *);
126285c0a3SHansol  Suh   PetscErrorCode (*regjac)(Tao, Vec, Mat, Mat, void *);
136285c0a3SHansol  Suh };
146285c0a3SHansol  Suh 
156285c0a3SHansol  Suh typedef struct {
166285c0a3SHansol  Suh   PETSCHEADER(struct _TaoADMMOps);
176285c0a3SHansol  Suh   Tao                    subsolverX, subsolverZ, parent;
186285c0a3SHansol  Suh   Vec                    residual, y, yold, y0, yhat, yhatold, constraint;
196285c0a3SHansol  Suh   Vec                    z, zold, Ax, Bz, Axold, Bzold, Bz0;
206285c0a3SHansol  Suh   Vec                    workLeft, workJacobianRight, workJacobianRight2; /*Ax,Bz,y,constraint are workJacobianRight sized. workLeft is solution sized */
216285c0a3SHansol  Suh   Mat                    Hx, Hxpre, Hz, Hzpre, ATA, BTB, JA, JApre, JB, JBpre;
226285c0a3SHansol  Suh   void                  *regobjgradP;
236285c0a3SHansol  Suh   void                  *reghessP;
246285c0a3SHansol  Suh   void                  *regjacobianP;
256285c0a3SHansol  Suh   void                  *misfitobjgradP;
266285c0a3SHansol  Suh   void                  *misfithessP;
276285c0a3SHansol  Suh   void                  *misfitjacobianP;
286285c0a3SHansol  Suh   PetscReal              gamma, last_misfit_val, last_reg_val, l1epsilon;
296285c0a3SHansol  Suh   PetscReal              lambda, mu, muold, orthval, mueps, tol, mumin;
306285c0a3SHansol  Suh   PetscReal              Bzdiffnorm, dualres, resnorm, const_norm;
316285c0a3SHansol  Suh   PetscReal              gatol_admm, catol_admm;
327f5c9be9SBarry Smith   PetscInt               T;                  /* adaptive iteration cutoff */
336285c0a3SHansol  Suh   PetscBool              xJI, zJI;           /* Bool to check whether A,B Jacobians are NULL-set identity */
346285c0a3SHansol  Suh   PetscBool              Hxchange, Hzchange; /* Bool to check whether Hx,Hz change wrt to x and z */
356285c0a3SHansol  Suh   PetscBool              Hxbool, Hzbool;     /* Bool to make sure Hessian gets updated only once for Hchange False case */
366285c0a3SHansol  Suh   TaoADMMUpdateType      update;             /* update policy for mu */
376285c0a3SHansol  Suh   TaoADMMRegularizerType regswitch;          /* regularization policy */
386285c0a3SHansol  Suh } TAO_ADMM;
39