xref: /petsc/src/tao/constrained/impls/ipm/pdipm.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
1*a4963045SJacob Faibussowitsch #pragma once
2aad13602SShrirang Abhyankar #include <petsc/private/taoimpl.h>
3aad13602SShrirang Abhyankar 
4aad13602SShrirang Abhyankar /*
5aad13602SShrirang Abhyankar  Context for Primal-Dual Interior-Point Method
6aad13602SShrirang Abhyankar  See the document pdipm.pdf
7aad13602SShrirang Abhyankar */
8aad13602SShrirang Abhyankar 
9aad13602SShrirang Abhyankar typedef struct {
10aad13602SShrirang Abhyankar   /* Sizes (n = local, N = global) */
11aad13602SShrirang Abhyankar   PetscInt nx, Nx;           /* Decision variables nx = nxfixed + nxub + nxlb + nxbox + nxfree */
12aad13602SShrirang Abhyankar   PetscInt nxfixed, Nxfixed; /* Fixed decision variables */
13aad13602SShrirang Abhyankar   PetscInt nxlb, Nxlb;       /* Decision variables with lower bounds only */
14aad13602SShrirang Abhyankar   PetscInt nxub, Nxub;       /* Decision variables with upper bounds only */
15aad13602SShrirang Abhyankar   PetscInt nxbox, Nxbox;     /* Decision variables with box constraints */
16aad13602SShrirang Abhyankar   PetscInt nxfree, Nxfree;   /* Free variables */
17aad13602SShrirang Abhyankar   PetscInt ng, Ng;           /* user equality constraints g(x) = 0. */
18aad13602SShrirang Abhyankar   PetscInt nh, Nh;           /* user inequality constraints h(x) >= 0. */
19aad13602SShrirang Abhyankar   PetscInt nce, Nce;         /* total equality constraints. nce = ng + nxfixed */
20aad13602SShrirang Abhyankar   PetscInt nci, Nci;         /* total inequality constraints nci = nh + nxlb + nxub + 2*nxbox */
21aad13602SShrirang Abhyankar   PetscInt n, N;             /* Big KKT system size n = nx + nce + 2*nci */
22aad13602SShrirang Abhyankar 
23aad13602SShrirang Abhyankar   /* Vectors */
24aad13602SShrirang Abhyankar   Vec X;       /* R^n   - Big KKT system vector [x; lambdae; lambdai; z] */
25aad13602SShrirang Abhyankar   Vec x;       /* R^nx - work vector, same layout as tao->solution */
26aad13602SShrirang Abhyankar   Vec lambdae; /* R^nce - vector, shares local arrays with X */
27aad13602SShrirang Abhyankar   Vec lambdai; /* R^nci - vector, shares local arrays with X */
28aad13602SShrirang Abhyankar   Vec z;       /* R^nci - vector, shares local arrays with X */
29aad13602SShrirang Abhyankar 
30aad13602SShrirang Abhyankar   /* Work vectors */
316aad120cSJose E. Roman   Vec lambdae_xfixed; /* Equality constraints lagrangian multiplier vector for fixed variables */
326aad120cSJose E. Roman   Vec lambdai_xb;     /* User inequality constraints lagrangian multiplier vector */
33aad13602SShrirang Abhyankar 
34aad13602SShrirang Abhyankar   /* Lagrangian equality and inequality Vec */
35aad13602SShrirang Abhyankar   Vec ce, ci; /* equality and inequality constraints */
36aad13602SShrirang Abhyankar 
37aad13602SShrirang Abhyankar   /* Offsets for subvectors */
38aad13602SShrirang Abhyankar   PetscInt off_lambdae, off_lambdai, off_z;
39aad13602SShrirang Abhyankar 
40aad13602SShrirang Abhyankar   /* Scalars */
41aad13602SShrirang Abhyankar   PetscReal L;     /* Lagrangian = f(x) - lambdae^T*ce(x) - lambdai^T*(ci(x) - z) - mu*sum_{i=1}^{Nci}(log(z_i)) */
42aad13602SShrirang Abhyankar   PetscReal gradL; /* gradient of L w.r.t. x */
43aad13602SShrirang Abhyankar 
44aad13602SShrirang Abhyankar   /* Matrices */
45aad13602SShrirang Abhyankar   Mat Jce_xfixed; /* Jacobian of equality constraints cebound(x) = J(nxfixed) */
46aad13602SShrirang Abhyankar   Mat Jci_xb;     /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */
47aad13602SShrirang Abhyankar   Mat K;          /* KKT matrix */
48aad13602SShrirang Abhyankar 
49aad13602SShrirang Abhyankar   /* Parameters */
50aad13602SShrirang Abhyankar   PetscReal mu;               /* Barrier parameter */
51aad13602SShrirang Abhyankar   PetscReal mu_update_factor; /* Multiplier for mu update */
5209ee8bb0SRylee Sundermann   PetscReal deltaw;
5309ee8bb0SRylee Sundermann   PetscReal lastdeltaw;
5409ee8bb0SRylee Sundermann   PetscReal deltac;
55aad13602SShrirang Abhyankar 
56aad13602SShrirang Abhyankar   /* Tolerances */
57aad13602SShrirang Abhyankar 
58aad13602SShrirang Abhyankar   /* Index sets for types of bounds on variables */
59aad13602SShrirang Abhyankar   IS isxub;    /* Finite upper bound only -inf < x < ub   */
60aad13602SShrirang Abhyankar   IS isxlb;    /* Finite lower bound only  lb <= x < inf  */
61aad13602SShrirang Abhyankar   IS isxfixed; /* Fixed variables          lb  = x = ub   */
62aad13602SShrirang Abhyankar   IS isxbox;   /* Boxed variables          lb <= x <= ub  */
63aad13602SShrirang Abhyankar   IS isxfree;  /* Free variables         -inf <= x <= inf */
64aad13602SShrirang Abhyankar 
65aad13602SShrirang Abhyankar   /* Index sets for PC fieldsplit */
66aad13602SShrirang Abhyankar   IS is1, is2;
67aad13602SShrirang Abhyankar 
68aad13602SShrirang Abhyankar   /* Options */
69aad13602SShrirang Abhyankar   PetscBool monitorkkt;          /* Monitor KKT */
70aad13602SShrirang Abhyankar   PetscReal push_init_slack;     /* Push initial slack variables (z) away from bounds */
71aad13602SShrirang Abhyankar   PetscReal push_init_lambdai;   /* Push initial inequality variables (lambdai) away from bounds */
72aad13602SShrirang Abhyankar   PetscBool solve_reduced_kkt;   /* Solve Reduced KKT with fieldsplit */
7312d688e0SRylee Sundermann   PetscBool solve_symmetric_kkt; /* Solve non-reduced symmetric KKT system */
7409ee8bb0SRylee Sundermann   PetscBool kkt_pd;              /* Add deltaw and deltac shifts to make KKT matrix positive definite */
75aad13602SShrirang Abhyankar 
76aad13602SShrirang Abhyankar   SNES snes;                                     /* Nonlinear solver */
77aad13602SShrirang Abhyankar   Mat  jac_equality_trans, jac_inequality_trans; /* working matrices */
78aad13602SShrirang Abhyankar 
79aad13602SShrirang Abhyankar   PetscReal obj; /* Objective function */
80aad13602SShrirang Abhyankar 
81aad13602SShrirang Abhyankar   /* Offsets for parallel assembly */
82aad13602SShrirang Abhyankar   PetscInt *nce_all;
83aad13602SShrirang Abhyankar } TAO_PDIPM;
84