xref: /petsc/src/tao/bound/impls/bnk/bnk.h (revision 080d2917cae42b0a43dd96fa397634a71dcc282e)
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>
9eb910715SAlp Dener 
10eb910715SAlp Dener typedef struct {
11eb910715SAlp Dener   Mat M;
12eb910715SAlp Dener   Vec W;
13eb910715SAlp Dener 
14eb910715SAlp Dener   Vec Xold;
15eb910715SAlp Dener   Vec Gold;
16eb910715SAlp Dener   Vec Diag;
17eb910715SAlp Dener 
18*080d2917SAlp Dener   /* Scalar values for the solution and step */
19*080d2917SAlp Dener   PetscReal fold, f, gnorm, dnorm;
20eb910715SAlp Dener 
21eb910715SAlp Dener   /* Parameters when updating the perturbation added to the Hessian matrix
22eb910715SAlp Dener      according to the following scheme:
23eb910715SAlp Dener 
24eb910715SAlp Dener      pert = sval;
25eb910715SAlp Dener 
26eb910715SAlp Dener      do until convergence
27eb910715SAlp Dener        shift Hessian by pert
28eb910715SAlp Dener        solve Newton system
29eb910715SAlp Dener 
30eb910715SAlp Dener        if (linear solver failed or did not compute a descent direction)
31eb910715SAlp Dener          use steepest descent direction and increase perturbation
32eb910715SAlp Dener 
33eb910715SAlp Dener          if (0 == pert)
34eb910715SAlp Dener            initialize perturbation
35eb910715SAlp Dener            pert = min(imax, max(imin, imfac * norm(G)))
36eb910715SAlp Dener          else
37eb910715SAlp Dener            increase perturbation
38eb910715SAlp Dener            pert = min(pmax, max(pgfac * pert, pmgfac * norm(G)))
39eb910715SAlp Dener          fi
40eb910715SAlp Dener        else
41eb910715SAlp Dener          use linear solver direction and decrease perturbation
42eb910715SAlp Dener 
43eb910715SAlp Dener          pert = min(psfac * pert, pmsfac * norm(G))
44eb910715SAlp Dener          if (pert < pmin)
45eb910715SAlp Dener            pert = 0
46eb910715SAlp Dener          fi
47eb910715SAlp Dener        fi
48eb910715SAlp Dener 
49eb910715SAlp Dener        perform line search
50eb910715SAlp Dener        function and gradient evaluation
51eb910715SAlp Dener        check convergence
52eb910715SAlp Dener      od
53eb910715SAlp Dener   */
54eb910715SAlp Dener   PetscReal sval;               /*  Starting perturbation value, default zero */
55eb910715SAlp Dener 
56eb910715SAlp Dener   PetscReal imin;               /*  Minimum perturbation added during initialization  */
57eb910715SAlp Dener   PetscReal imax;               /*  Maximum perturbation added during initialization */
58eb910715SAlp Dener   PetscReal imfac;              /*  Merit function factor during initialization */
59eb910715SAlp Dener 
60eb910715SAlp Dener   PetscReal pert;               /*  Current perturbation value */
61eb910715SAlp Dener   PetscReal pmin;               /*  Minimim perturbation value */
62eb910715SAlp Dener   PetscReal pmax;               /*  Maximum perturbation value */
63eb910715SAlp Dener   PetscReal pgfac;              /*  Perturbation growth factor */
64eb910715SAlp Dener   PetscReal psfac;              /*  Perturbation shrink factor */
65eb910715SAlp Dener   PetscReal pmgfac;             /*  Merit function growth factor */
66eb910715SAlp Dener   PetscReal pmsfac;             /*  Merit function shrink factor */
67eb910715SAlp Dener 
68eb910715SAlp Dener   /* Parameters when updating the trust-region radius based on steplength
69eb910715SAlp Dener      if   step < nu1            (very bad step)
70eb910715SAlp Dener        radius = omega1 * min(norm(d), radius)
71eb910715SAlp Dener      elif step < nu2            (bad step)
72eb910715SAlp Dener        radius = omega2 * min(norm(d), radius)
73eb910715SAlp Dener      elif step < nu3            (okay step)
74eb910715SAlp Dener        radius = omega3 * radius;
75eb910715SAlp Dener      elif step < nu4            (good step)
76eb910715SAlp Dener        radius = max(omega4 * norm(d), radius)
77eb910715SAlp Dener      else                       (very good step)
78eb910715SAlp Dener        radius = max(omega5 * norm(d), radius)
79eb910715SAlp Dener      fi
80eb910715SAlp Dener   */
81eb910715SAlp Dener   PetscReal nu1;                /*  used to compute trust-region radius */
82eb910715SAlp Dener   PetscReal nu2;                /*  used to compute trust-region radius */
83eb910715SAlp Dener   PetscReal nu3;                /*  used to compute trust-region radius */
84eb910715SAlp Dener   PetscReal nu4;                /*  used to compute trust-region radius */
85eb910715SAlp Dener 
86eb910715SAlp Dener   PetscReal omega1;             /*  factor used for trust-region update */
87eb910715SAlp Dener   PetscReal omega2;             /*  factor used for trust-region update */
88eb910715SAlp Dener   PetscReal omega3;             /*  factor used for trust-region update */
89eb910715SAlp Dener   PetscReal omega4;             /*  factor used for trust-region update */
90eb910715SAlp Dener   PetscReal omega5;             /*  factor used for trust-region update */
91eb910715SAlp Dener 
92eb910715SAlp Dener   /* Parameters when updating the trust-region radius based on reduction
93eb910715SAlp Dener 
94eb910715SAlp Dener      kappa = ared / pred
95eb910715SAlp Dener      if   kappa < eta1          (very bad step)
96eb910715SAlp Dener        radius = alpha1 * min(norm(d), radius)
97eb910715SAlp Dener      elif kappa < eta2          (bad step)
98eb910715SAlp Dener        radius = alpha2 * min(norm(d), radius)
99eb910715SAlp Dener      elif kappa < eta3          (okay step)
100eb910715SAlp Dener        radius = alpha3 * radius;
101eb910715SAlp Dener      elif kappa < eta4          (good step)
102eb910715SAlp Dener        radius = max(alpha4 * norm(d), radius)
103eb910715SAlp Dener      else                       (very good step)
104eb910715SAlp Dener        radius = max(alpha5 * norm(d), radius)
105eb910715SAlp Dener      fi
106eb910715SAlp Dener   */
107eb910715SAlp Dener   PetscReal eta1;               /*  used to compute trust-region radius */
108eb910715SAlp Dener   PetscReal eta2;               /*  used to compute trust-region radius */
109eb910715SAlp Dener   PetscReal eta3;               /*  used to compute trust-region radius */
110eb910715SAlp Dener   PetscReal eta4;               /*  used to compute trust-region radius */
111eb910715SAlp Dener 
112eb910715SAlp Dener   PetscReal alpha1;             /*  factor used for trust-region update */
113eb910715SAlp Dener   PetscReal alpha2;             /*  factor used for trust-region update */
114eb910715SAlp Dener   PetscReal alpha3;             /*  factor used for trust-region update */
115eb910715SAlp Dener   PetscReal alpha4;             /*  factor used for trust-region update */
116eb910715SAlp Dener   PetscReal alpha5;             /*  factor used for trust-region update */
117eb910715SAlp Dener 
118eb910715SAlp Dener   /* Parameters when updating the trust-region radius based on interpolation
119eb910715SAlp Dener 
120eb910715SAlp Dener      kappa = ared / pred
121eb910715SAlp Dener      if   kappa >= 1.0 - mu1    (very good step)
122eb910715SAlp Dener        choose tau in [gamma3, gamma4]
123eb910715SAlp Dener        radius = max(tau * norm(d), radius)
124eb910715SAlp Dener      elif kappa >= 1.0 - mu2    (good step)
125eb910715SAlp Dener        choose tau in [gamma2, gamma3]
126eb910715SAlp Dener        if (tau >= 1.0)
127eb910715SAlp Dener          radius = max(tau * norm(d), radius)
128eb910715SAlp Dener        else
129eb910715SAlp Dener          radius = tau * min(norm(d), radius)
130eb910715SAlp Dener        fi
131eb910715SAlp Dener      else                       (bad step)
132eb910715SAlp Dener        choose tau in [gamma1, 1.0]
133eb910715SAlp Dener        radius = tau * min(norm(d), radius)
134eb910715SAlp Dener      fi
135eb910715SAlp Dener   */
136eb910715SAlp Dener   PetscReal mu1;                /*  used for model agreement in interpolation */
137eb910715SAlp Dener   PetscReal mu2;                /*  used for model agreement in interpolation */
138eb910715SAlp Dener 
139eb910715SAlp Dener   PetscReal gamma1;             /*  factor used for interpolation */
140eb910715SAlp Dener   PetscReal gamma2;             /*  factor used for interpolation */
141eb910715SAlp Dener   PetscReal gamma3;             /*  factor used for interpolation */
142eb910715SAlp Dener   PetscReal gamma4;             /*  factor used for interpolation */
143eb910715SAlp Dener 
144eb910715SAlp Dener   PetscReal theta;              /*  factor used for interpolation */
145eb910715SAlp Dener 
146eb910715SAlp Dener   /*  Parameters when initializing trust-region radius based on interpolation */
147eb910715SAlp Dener   PetscReal mu1_i;              /*  used for model agreement in interpolation */
148eb910715SAlp Dener   PetscReal mu2_i;              /*  used for model agreement in interpolation */
149eb910715SAlp Dener 
150eb910715SAlp Dener   PetscReal gamma1_i;           /*  factor used for interpolation */
151eb910715SAlp Dener   PetscReal gamma2_i;           /*  factor used for interpolation */
152eb910715SAlp Dener   PetscReal gamma3_i;           /*  factor used for interpolation */
153eb910715SAlp Dener   PetscReal gamma4_i;           /*  factor used for interpolation */
154eb910715SAlp Dener 
155eb910715SAlp Dener   PetscReal theta_i;            /*  factor used for interpolation */
156eb910715SAlp Dener 
157eb910715SAlp Dener   /*  Other parameters */
158eb910715SAlp Dener   PetscReal min_radius;         /*  lower bound on initial radius value */
159eb910715SAlp Dener   PetscReal max_radius;         /*  upper bound on trust region radius */
160eb910715SAlp Dener   PetscReal epsilon;            /*  tolerance used when computing ared/pred */
161eb910715SAlp Dener 
162eb910715SAlp Dener   PetscInt newt;                /*  Newton directions attempted */
163eb910715SAlp Dener   PetscInt bfgs;                /*  BFGS directions attempted */
164eb910715SAlp Dener   PetscInt sgrad;               /*  Scaled gradient directions attempted */
165eb910715SAlp Dener   PetscInt grad;                /*  Gradient directions attempted */
166eb910715SAlp Dener 
167eb910715SAlp Dener   PetscInt pc_type;             /*  Preconditioner for the code */
168eb910715SAlp Dener   PetscInt bfgs_scale_type;     /*  Scaling matrix to used for the bfgs preconditioner */
169eb910715SAlp Dener   PetscInt init_type;           /*  Trust-region initialization method */
170eb910715SAlp Dener   PetscInt update_type;         /*  Trust-region update method */
171eb910715SAlp Dener 
172eb910715SAlp Dener   PetscInt ksp_atol;
173eb910715SAlp Dener   PetscInt ksp_rtol;
174eb910715SAlp Dener   PetscInt ksp_ctol;
175eb910715SAlp Dener   PetscInt ksp_negc;
176eb910715SAlp Dener   PetscInt ksp_dtol;
177eb910715SAlp Dener   PetscInt ksp_iter;
178eb910715SAlp Dener   PetscInt ksp_othr;
179eb910715SAlp Dener   PetscBool is_nash, is_stcg, is_gltr;
180eb910715SAlp Dener } TAO_BNK;
181eb910715SAlp Dener 
182eb910715SAlp Dener #endif /* if !defined(__TAO_BNK_H) */
183eb910715SAlp Dener 
184eb910715SAlp Dener #define BNK_NEWTON              0
185eb910715SAlp Dener #define BNK_BFGS                1
186eb910715SAlp Dener #define BNK_SCALED_GRADIENT     2
187eb910715SAlp Dener #define BNK_GRADIENT            3
188eb910715SAlp Dener 
189eb910715SAlp Dener #define BNK_PC_NONE     0
190eb910715SAlp Dener #define BNK_PC_AHESS    1
191eb910715SAlp Dener #define BNK_PC_BFGS     2
192eb910715SAlp Dener #define BNK_PC_PETSC    3
193eb910715SAlp Dener #define BNK_PC_TYPES    4
194eb910715SAlp Dener 
195eb910715SAlp Dener #define BFGS_SCALE_AHESS        0
196eb910715SAlp Dener #define BFGS_SCALE_PHESS        1
197eb910715SAlp Dener #define BFGS_SCALE_BFGS         2
198eb910715SAlp Dener #define BFGS_SCALE_TYPES        3
199eb910715SAlp Dener 
200eb910715SAlp Dener #define BNK_INIT_CONSTANT         0
201eb910715SAlp Dener #define BNK_INIT_DIRECTION        1
202eb910715SAlp Dener #define BNK_INIT_INTERPOLATION    2
203eb910715SAlp Dener #define BNK_INIT_TYPES            3
204eb910715SAlp Dener 
205eb910715SAlp Dener #define BNK_UPDATE_STEP           0
206eb910715SAlp Dener #define BNK_UPDATE_REDUCTION      1
207eb910715SAlp Dener #define BNK_UPDATE_INTERPOLATION  2
208eb910715SAlp Dener #define BNK_UPDATE_TYPES          3
209eb910715SAlp Dener 
210eb910715SAlp Dener static const char *BNK_PC[64] = {"none", "ahess", "bfgs", "petsc"};
211eb910715SAlp Dener 
212eb910715SAlp Dener static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
213eb910715SAlp Dener 
214eb910715SAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
215eb910715SAlp Dener 
216eb910715SAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
217eb910715SAlp Dener 
218eb910715SAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BNK(Tao);
219eb910715SAlp Dener 
220eb910715SAlp Dener PETSC_INTERN PetscErrorCode MatLMVMSolveShell(PC, Vec, Vec);
221eb910715SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKInitialize(Tao);
222eb910715SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKComputeStep(Tao, PetscInt*);
223*080d2917SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKUpdateTrustRadius(Tao, PetscReal, PetscReal, PetscInt, PetscBool*);
224