xref: /petsc/src/tao/bound/impls/bnk/bnk.h (revision 0875260307e97403c513ee2ad2ca062f01ff5a8a)
1eb910715SAlp Dener /*
2eb910715SAlp Dener Context for bounded Newton-Krylov type optimization algorithms
3eb910715SAlp Dener */
4eb910715SAlp Dener 
5eb910715SAlp Dener #if !defined(__TAO_BNK_H)
6eb910715SAlp Dener #define __TAO_BNK_H
7eb910715SAlp Dener #include <petsc/private/taoimpl.h>
8eb910715SAlp Dener #include <../src/tao/matrix/lmvmmat.h>
9c0f10754SAlp Dener #include <../src/tao/bound/impls/bncg/bncg.h>
10eb910715SAlp Dener 
11eb910715SAlp Dener typedef struct {
12c0f10754SAlp Dener   /* Embedded TAOBNCG */
13c0f10754SAlp Dener   Tao bncg;
14c0f10754SAlp Dener   TAO_BNCG *bncg_ctx;
15e031d6f5SAlp Dener   PetscInt max_cg_its, tot_cg_its;
16c0f10754SAlp Dener   Vec bncg_sol;
17c0f10754SAlp Dener 
18c0f10754SAlp Dener   /* Allocated vectors */
19c0f10754SAlp Dener   Vec W, Xwork, Gwork, Xold, Gold;
2009164190SAlp Dener   Vec unprojected_gradient, unprojected_gradient_old;
21c0f10754SAlp Dener   Vec Diag_min, Diag_max;
22c0f10754SAlp Dener 
23c0f10754SAlp Dener   /* Unallocated matrices and vectors */
24c0f10754SAlp Dener   Mat H_inactive, Hpre_inactive, M;
25c0f10754SAlp Dener   Vec Diag, Diag_red, X_inactive, G_inactive, inactive_work, active_work;
262f75a4aaSAlp Dener   IS  inactive_idx, active_idx, active_lower, active_upper, active_fixed;
27eb910715SAlp Dener 
28080d2917SAlp Dener   /* Scalar values for the solution and step */
29080d2917SAlp Dener   PetscReal fold, f, gnorm, dnorm;
30eb910715SAlp Dener 
312f75a4aaSAlp Dener   /* Parameters for active set estimation */
320a4511e9SAlp Dener   PetscReal as_tol;
330a4511e9SAlp Dener   PetscReal as_step;
342f75a4aaSAlp Dener 
35eb910715SAlp Dener   /* Parameters when updating the perturbation added to the Hessian matrix
36eb910715SAlp Dener      according to the following scheme:
37eb910715SAlp Dener 
38eb910715SAlp Dener      pert = sval;
39eb910715SAlp Dener 
40eb910715SAlp Dener      do until convergence
41eb910715SAlp Dener        shift Hessian by pert
42eb910715SAlp Dener        solve Newton system
43eb910715SAlp Dener 
44eb910715SAlp Dener        if (linear solver failed or did not compute a descent direction)
45eb910715SAlp Dener          use steepest descent direction and increase perturbation
46eb910715SAlp Dener 
47eb910715SAlp Dener          if (0 == pert)
48eb910715SAlp Dener            initialize perturbation
49eb910715SAlp Dener            pert = min(imax, max(imin, imfac * norm(G)))
50eb910715SAlp Dener          else
51eb910715SAlp Dener            increase perturbation
52eb910715SAlp Dener            pert = min(pmax, max(pgfac * pert, pmgfac * norm(G)))
53eb910715SAlp Dener          fi
54eb910715SAlp Dener        else
55eb910715SAlp Dener          use linear solver direction and decrease perturbation
56eb910715SAlp Dener 
57eb910715SAlp Dener          pert = min(psfac * pert, pmsfac * norm(G))
58eb910715SAlp Dener          if (pert < pmin)
59eb910715SAlp Dener            pert = 0
60eb910715SAlp Dener          fi
61eb910715SAlp Dener        fi
62eb910715SAlp Dener 
63eb910715SAlp Dener        perform line search
64eb910715SAlp Dener        function and gradient evaluation
65eb910715SAlp Dener        check convergence
66eb910715SAlp Dener      od
67eb910715SAlp Dener   */
68eb910715SAlp Dener   PetscReal sval;               /*  Starting perturbation value, default zero */
69eb910715SAlp Dener 
70eb910715SAlp Dener   PetscReal imin;               /*  Minimum perturbation added during initialization  */
71eb910715SAlp Dener   PetscReal imax;               /*  Maximum perturbation added during initialization */
72eb910715SAlp Dener   PetscReal imfac;              /*  Merit function factor during initialization */
73eb910715SAlp Dener 
74eb910715SAlp Dener   PetscReal pert;               /*  Current perturbation value */
75eb910715SAlp Dener   PetscReal pmin;               /*  Minimim perturbation value */
76eb910715SAlp Dener   PetscReal pmax;               /*  Maximum perturbation value */
77eb910715SAlp Dener   PetscReal pgfac;              /*  Perturbation growth factor */
78eb910715SAlp Dener   PetscReal psfac;              /*  Perturbation shrink factor */
79eb910715SAlp Dener   PetscReal pmgfac;             /*  Merit function growth factor */
80eb910715SAlp Dener   PetscReal pmsfac;             /*  Merit function shrink factor */
81eb910715SAlp Dener 
82eb910715SAlp Dener   /* Parameters when updating the trust-region radius based on steplength
83eb910715SAlp Dener      if   step < nu1            (very bad step)
84eb910715SAlp Dener        radius = omega1 * min(norm(d), radius)
85eb910715SAlp Dener      elif step < nu2            (bad step)
86eb910715SAlp Dener        radius = omega2 * min(norm(d), radius)
87eb910715SAlp Dener      elif step < nu3            (okay step)
88eb910715SAlp Dener        radius = omega3 * radius;
89eb910715SAlp Dener      elif step < nu4            (good step)
90eb910715SAlp Dener        radius = max(omega4 * norm(d), radius)
91eb910715SAlp Dener      else                       (very good step)
92eb910715SAlp Dener        radius = max(omega5 * norm(d), radius)
93eb910715SAlp Dener      fi
94eb910715SAlp Dener   */
95eb910715SAlp Dener   PetscReal nu1;                /*  used to compute trust-region radius */
96eb910715SAlp Dener   PetscReal nu2;                /*  used to compute trust-region radius */
97eb910715SAlp Dener   PetscReal nu3;                /*  used to compute trust-region radius */
98eb910715SAlp Dener   PetscReal nu4;                /*  used to compute trust-region radius */
99eb910715SAlp Dener 
100eb910715SAlp Dener   PetscReal omega1;             /*  factor used for trust-region update */
101eb910715SAlp Dener   PetscReal omega2;             /*  factor used for trust-region update */
102eb910715SAlp Dener   PetscReal omega3;             /*  factor used for trust-region update */
103eb910715SAlp Dener   PetscReal omega4;             /*  factor used for trust-region update */
104eb910715SAlp Dener   PetscReal omega5;             /*  factor used for trust-region update */
105eb910715SAlp Dener 
106eb910715SAlp Dener   /* Parameters when updating the trust-region radius based on reduction
107eb910715SAlp Dener 
108eb910715SAlp Dener      kappa = ared / pred
109eb910715SAlp Dener      if   kappa < eta1          (very bad step)
110eb910715SAlp Dener        radius = alpha1 * min(norm(d), radius)
111eb910715SAlp Dener      elif kappa < eta2          (bad step)
112eb910715SAlp Dener        radius = alpha2 * min(norm(d), radius)
113eb910715SAlp Dener      elif kappa < eta3          (okay step)
114eb910715SAlp Dener        radius = alpha3 * radius;
115eb910715SAlp Dener      elif kappa < eta4          (good step)
116eb910715SAlp Dener        radius = max(alpha4 * norm(d), radius)
117eb910715SAlp Dener      else                       (very good step)
118eb910715SAlp Dener        radius = max(alpha5 * norm(d), radius)
119eb910715SAlp Dener      fi
120eb910715SAlp Dener   */
121eb910715SAlp Dener   PetscReal eta1;               /*  used to compute trust-region radius */
122eb910715SAlp Dener   PetscReal eta2;               /*  used to compute trust-region radius */
123eb910715SAlp Dener   PetscReal eta3;               /*  used to compute trust-region radius */
124eb910715SAlp Dener   PetscReal eta4;               /*  used to compute trust-region radius */
125eb910715SAlp Dener 
126eb910715SAlp Dener   PetscReal alpha1;             /*  factor used for trust-region update */
127eb910715SAlp Dener   PetscReal alpha2;             /*  factor used for trust-region update */
128eb910715SAlp Dener   PetscReal alpha3;             /*  factor used for trust-region update */
129eb910715SAlp Dener   PetscReal alpha4;             /*  factor used for trust-region update */
130eb910715SAlp Dener   PetscReal alpha5;             /*  factor used for trust-region update */
131eb910715SAlp Dener 
132eb910715SAlp Dener   /* Parameters when updating the trust-region radius based on interpolation
133eb910715SAlp Dener 
134eb910715SAlp Dener      kappa = ared / pred
135eb910715SAlp Dener      if   kappa >= 1.0 - mu1    (very good step)
136eb910715SAlp Dener        choose tau in [gamma3, gamma4]
137eb910715SAlp Dener        radius = max(tau * norm(d), radius)
138eb910715SAlp Dener      elif kappa >= 1.0 - mu2    (good step)
139eb910715SAlp Dener        choose tau in [gamma2, gamma3]
140eb910715SAlp Dener        if (tau >= 1.0)
141eb910715SAlp Dener          radius = max(tau * norm(d), radius)
142eb910715SAlp Dener        else
143eb910715SAlp Dener          radius = tau * min(norm(d), radius)
144eb910715SAlp Dener        fi
145eb910715SAlp Dener      else                       (bad step)
146eb910715SAlp Dener        choose tau in [gamma1, 1.0]
147eb910715SAlp Dener        radius = tau * min(norm(d), radius)
148eb910715SAlp Dener      fi
149eb910715SAlp Dener   */
150eb910715SAlp Dener   PetscReal mu1;                /*  used for model agreement in interpolation */
151eb910715SAlp Dener   PetscReal mu2;                /*  used for model agreement in interpolation */
152eb910715SAlp Dener 
153eb910715SAlp Dener   PetscReal gamma1;             /*  factor used for interpolation */
154eb910715SAlp Dener   PetscReal gamma2;             /*  factor used for interpolation */
155eb910715SAlp Dener   PetscReal gamma3;             /*  factor used for interpolation */
156eb910715SAlp Dener   PetscReal gamma4;             /*  factor used for interpolation */
157eb910715SAlp Dener 
158eb910715SAlp Dener   PetscReal theta;              /*  factor used for interpolation */
159eb910715SAlp Dener 
160eb910715SAlp Dener   /*  Parameters when initializing trust-region radius based on interpolation */
161eb910715SAlp Dener   PetscReal mu1_i;              /*  used for model agreement in interpolation */
162eb910715SAlp Dener   PetscReal mu2_i;              /*  used for model agreement in interpolation */
163eb910715SAlp Dener 
164eb910715SAlp Dener   PetscReal gamma1_i;           /*  factor used for interpolation */
165eb910715SAlp Dener   PetscReal gamma2_i;           /*  factor used for interpolation */
166eb910715SAlp Dener   PetscReal gamma3_i;           /*  factor used for interpolation */
167eb910715SAlp Dener   PetscReal gamma4_i;           /*  factor used for interpolation */
168eb910715SAlp Dener 
169eb910715SAlp Dener   PetscReal theta_i;            /*  factor used for interpolation */
170eb910715SAlp Dener 
171eb910715SAlp Dener   /*  Other parameters */
172eb910715SAlp Dener   PetscReal min_radius;         /*  lower bound on initial radius value */
173eb910715SAlp Dener   PetscReal max_radius;         /*  upper bound on trust region radius */
174eb910715SAlp Dener   PetscReal epsilon;            /*  tolerance used when computing ared/pred */
17562675beeSAlp Dener   PetscReal dmin, dmax;         /*  upper and lower bounds for the Hessian diagonal vector */
176eb910715SAlp Dener 
177eb910715SAlp Dener   PetscInt newt;                /*  Newton directions attempted */
178eb910715SAlp Dener   PetscInt bfgs;                /*  BFGS directions attempted */
179eb910715SAlp Dener   PetscInt sgrad;               /*  Scaled gradient directions attempted */
180eb910715SAlp Dener   PetscInt grad;                /*  Gradient directions attempted */
181eb910715SAlp Dener 
18262675beeSAlp Dener   PetscInt as_type;             /*   Active set estimation method */
183eb910715SAlp Dener   PetscInt pc_type;             /*  Preconditioner for the code */
184eb910715SAlp Dener   PetscInt bfgs_scale_type;     /*  Scaling matrix to used for the bfgs preconditioner */
185eb910715SAlp Dener   PetscInt init_type;           /*  Trust-region initialization method */
186eb910715SAlp Dener   PetscInt update_type;         /*  Trust-region update method */
187eb910715SAlp Dener 
1882f75a4aaSAlp Dener   /* Trackers for KSP solution type and convergence reasons */
189eb910715SAlp Dener   PetscInt ksp_atol;
190eb910715SAlp Dener   PetscInt ksp_rtol;
191eb910715SAlp Dener   PetscInt ksp_ctol;
192eb910715SAlp Dener   PetscInt ksp_negc;
193eb910715SAlp Dener   PetscInt ksp_dtol;
194eb910715SAlp Dener   PetscInt ksp_iter;
195eb910715SAlp Dener   PetscInt ksp_othr;
196eb910715SAlp Dener   PetscBool is_nash, is_stcg, is_gltr;
197eb910715SAlp Dener } TAO_BNK;
198eb910715SAlp Dener 
199eb910715SAlp Dener #endif /* if !defined(__TAO_BNK_H) */
200eb910715SAlp Dener 
201eb910715SAlp Dener #define BNK_NEWTON              0
202eb910715SAlp Dener #define BNK_BFGS                1
203eb910715SAlp Dener #define BNK_SCALED_GRADIENT     2
204eb910715SAlp Dener #define BNK_GRADIENT            3
205eb910715SAlp Dener 
206eb910715SAlp Dener #define BNK_PC_NONE     0
207eb910715SAlp Dener #define BNK_PC_AHESS    1
208eb910715SAlp Dener #define BNK_PC_BFGS     2
209eb910715SAlp Dener #define BNK_PC_PETSC    3
210eb910715SAlp Dener #define BNK_PC_TYPES    4
211eb910715SAlp Dener 
212eb910715SAlp Dener #define BFGS_SCALE_AHESS        0
213eb910715SAlp Dener #define BFGS_SCALE_PHESS        1
214eb910715SAlp Dener #define BFGS_SCALE_BFGS         2
215eb910715SAlp Dener #define BFGS_SCALE_TYPES        3
216eb910715SAlp Dener 
217eb910715SAlp Dener #define BNK_INIT_CONSTANT         0
218eb910715SAlp Dener #define BNK_INIT_DIRECTION        1
219eb910715SAlp Dener #define BNK_INIT_INTERPOLATION    2
220eb910715SAlp Dener #define BNK_INIT_TYPES            3
221eb910715SAlp Dener 
222eb910715SAlp Dener #define BNK_UPDATE_STEP           0
223eb910715SAlp Dener #define BNK_UPDATE_REDUCTION      1
224eb910715SAlp Dener #define BNK_UPDATE_INTERPOLATION  2
225eb910715SAlp Dener #define BNK_UPDATE_TYPES          3
226eb910715SAlp Dener 
2272f75a4aaSAlp Dener #define BNK_AS_NONE             0
2282f75a4aaSAlp Dener #define BNK_AS_BERTSEKAS        1
2292f75a4aaSAlp Dener #define BNK_AS_TYPES            2
2302f75a4aaSAlp Dener 
231eb910715SAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BNK(Tao);
2329b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoSetUp_BNK(Tao);
233eb910715SAlp Dener 
234eb910715SAlp Dener PETSC_INTERN PetscErrorCode MatLMVMSolveShell(PC, Vec, Vec);
235c0f10754SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKInitialize(Tao, PetscInt, PetscBool*);
236*08752603SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKEstimateActiveSet(Tao, PetscInt);
23762675beeSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKComputeHessian(Tao);
2382f75a4aaSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKBoundStep(Tao, Vec);
239c0f10754SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKTakeCGSteps(Tao, PetscBool*);
24062675beeSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKComputeStep(Tao, PetscBool, KSPConvergedReason*);
2415e9b73cbSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKRecomputePred(Tao, Vec, PetscReal*);
242e465cd6fSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKSafeguardStep(Tao, KSPConvergedReason, PetscInt*);
243c14b763aSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKPerformLineSearch(Tao, PetscInt, PetscReal*, TaoLineSearchConvergedReason*);
24428017e9fSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKUpdateTrustRadius(Tao, PetscReal, PetscReal, PetscInt, PetscInt, PetscBool*);
24562675beeSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKAddStepCounts(Tao, PetscInt);
246