xref: /petsc/src/tao/constrained/impls/almm/almm.h (revision d756bedd70a89ca052be956bccd75c5761cb2ab4)
1*a4963045SJacob Faibussowitsch #pragma once
2661095bbSAlp Dener #include <petsc/private/taoimpl.h>
3661095bbSAlp Dener 
4661095bbSAlp Dener typedef struct {
5661095bbSAlp Dener   Tao subsolver, parent;          /* subsolver for aug-lag subproblem */
6661095bbSAlp Dener   PetscErrorCode (*sub_obj)(Tao); /* subsolver objective function */
7661095bbSAlp Dener   TaoALMMType type;               /* subsolver objective type */
8661095bbSAlp Dener 
9661095bbSAlp Dener   IS         *Pis, *Yis;           /* index sets to separate primal and dual vector spaces */
10661095bbSAlp Dener   VecScatter *Pscatter, *Yscatter; /* scatter objects to write into combined vector spaces */
11661095bbSAlp Dener 
12661095bbSAlp Dener   Mat  Ae, Ai;                              /* aliased constraint Jacobians (do not destroy!) */
13661095bbSAlp Dener   Vec  Px, LgradX, Ce, Ci, G;               /* aliased vectors (do not destroy!) */
14661095bbSAlp Dener   Vec  Ps, LgradS, Yi, Ye;                  /* sub-vectors for primal variables */
15661095bbSAlp Dener   Vec *Parr, P, PL, PU, *Yarr, Y, C;        /* arrays and vectors for combined vector spaces */
169a09327dSAlp Dener   Vec  Psub, Xwork, Cework, Ciwork, Cizero; /* work vectors */
17661095bbSAlp Dener 
18661095bbSAlp Dener   PetscReal Lval, fval, gnorm, cnorm, cenorm, cinorm, cnorm_old; /* scalar variables */
19661095bbSAlp Dener   PetscReal mu0, mu, mu_fac, mu_pow_good, mu_pow_bad;            /* penalty parameters */
20661095bbSAlp Dener   PetscReal ytol0, ytol, gtol0, gtol;                            /* convergence parameters */
21661095bbSAlp Dener   PetscReal mu_max, ye_min, yi_min, ye_max, yi_max;              /* parameter safeguards */
22661095bbSAlp Dener } TAO_ALMM;
23661095bbSAlp Dener 
24661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMGetType_Private(Tao, TaoALMMType *);
25661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMSetType_Private(Tao, TaoALMMType);
26661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMGetSubsolver_Private(Tao, Tao *);
27661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMSetSubsolver_Private(Tao, Tao);
28661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMGetMultipliers_Private(Tao, Vec *);
29661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMSetMultipliers_Private(Tao, Vec);
30661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMGetPrimalIS_Private(Tao, IS *, IS *);
31661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMGetDualIS_Private(Tao, IS *, IS *);
327721a36fSStefano Zampini PETSC_INTERN PetscErrorCode TaoALMMSubsolverObjective_Private(Tao, Vec, PetscReal *, void *);
33661095bbSAlp Dener PETSC_INTERN PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao, Vec, PetscReal *, Vec, void *);
34