xref: /petsc/src/tao/unconstrained/impls/ntr/ntrimpl.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
1fb90e4d1STodd Munson /*
2fb90e4d1STodd Munson   Context for a Newton trust region method (unconstrained minimization)
3fb90e4d1STodd Munson */
4fb90e4d1STodd Munson 
5*a4963045SJacob Faibussowitsch #pragma once
6fb90e4d1STodd Munson #include <petsc/private/taoimpl.h>
7fb90e4d1STodd Munson 
8fb90e4d1STodd Munson typedef struct {
9fb90e4d1STodd Munson   Mat M;
100c51296cSAlp Dener   PC  bfgs_pre;
11fb90e4d1STodd Munson 
12fb90e4d1STodd Munson   Vec D;
13fb90e4d1STodd Munson   Vec W;
14fb90e4d1STodd Munson 
15fb90e4d1STodd Munson   PetscReal radius;
16fb90e4d1STodd Munson 
17fb90e4d1STodd Munson   /* Parameters when updating the trust-region radius based on reduction
18fb90e4d1STodd Munson 
19fb90e4d1STodd Munson      kappa = ared / pred
20fb90e4d1STodd Munson      if   kappa < eta1          (very bad step)
21fb90e4d1STodd Munson        radius = alpha1 * min(norm(d), radius)
22fb90e4d1STodd Munson      elif kappa < eta2          (bad step)
23fb90e4d1STodd Munson        radius = alpha2 * min(norm(d), radius)
24fb90e4d1STodd Munson      elif kappa < eta3          (okay step)
25fb90e4d1STodd Munson        radius = alpha3 * radius;
26fb90e4d1STodd Munson      elif kappa < eta4          (good step)
27fb90e4d1STodd Munson        radius = max(alpha4 * norm(d), radius)
28fb90e4d1STodd Munson      else                       (very good step)
29fb90e4d1STodd Munson        radius = max(alpha5 * norm(d), radius)
30fb90e4d1STodd Munson      fi
31fb90e4d1STodd Munson   */
32fb90e4d1STodd Munson 
33fb90e4d1STodd Munson   PetscReal eta1; /*  used to compute trust-region radius */
34fb90e4d1STodd Munson   PetscReal eta2; /*  used to compute trust-region radius */
35fb90e4d1STodd Munson   PetscReal eta3; /*  used to compute trust-region radius */
36fb90e4d1STodd Munson   PetscReal eta4; /*  used to compute trust-region radius */
37fb90e4d1STodd Munson 
38fb90e4d1STodd Munson   PetscReal alpha1; /*  factor used for trust-region update */
39fb90e4d1STodd Munson   PetscReal alpha2; /*  factor used for trust-region update */
40fb90e4d1STodd Munson   PetscReal alpha3; /*  factor used for trust-region update */
41fb90e4d1STodd Munson   PetscReal alpha4; /*  factor used for trust-region update */
42fb90e4d1STodd Munson   PetscReal alpha5; /*  factor used for trust-region update */
43fb90e4d1STodd Munson 
44fb90e4d1STodd Munson   /* Parameters when updating the trust-region radius based on interpolation
45fb90e4d1STodd Munson 
46fb90e4d1STodd Munson      kappa = ared / pred
47fb90e4d1STodd Munson      if   kappa >= 1.0 - mu1    (very good step)
48fb90e4d1STodd Munson        choose tau in [gamma3, gamma4]
49fb90e4d1STodd Munson        radius = max(tau * norm(d), radius)
50fb90e4d1STodd Munson      elif kappa >= 1.0 - mu2    (good step)
51fb90e4d1STodd Munson        choose tau in [gamma2, gamma3]
52fb90e4d1STodd Munson        if (tau >= 1.0)
53fb90e4d1STodd Munson          radius = max(tau * norm(d), radius)
54fb90e4d1STodd Munson        else
55fb90e4d1STodd Munson          radius = tau * min(norm(d), radius)
56fb90e4d1STodd Munson        fi
57fb90e4d1STodd Munson      else                       (bad step)
58fb90e4d1STodd Munson        choose tau in [gamma1, 1.0]
59fb90e4d1STodd Munson        radius = tau * min(norm(d), radius)
60fb90e4d1STodd Munson      fi
61fb90e4d1STodd Munson   */
62fb90e4d1STodd Munson 
63fb90e4d1STodd Munson   PetscReal mu1; /*  used for model agreement in radius update */
64fb90e4d1STodd Munson   PetscReal mu2; /*  used for model agreement in radius update */
65fb90e4d1STodd Munson 
66fb90e4d1STodd Munson   PetscReal gamma1; /*  factor used for radius update */
67fb90e4d1STodd Munson   PetscReal gamma2; /*  factor used for radius update */
68fb90e4d1STodd Munson   PetscReal gamma3; /*  factor used for radius update */
69fb90e4d1STodd Munson   PetscReal gamma4; /*  factor used for radius update */
70fb90e4d1STodd Munson 
71fb90e4d1STodd Munson   PetscReal theta; /*  factor used for radius update */
72fb90e4d1STodd Munson 
73fb90e4d1STodd Munson   /* Parameters when initializing trust-region radius based on interpolation */
74fb90e4d1STodd Munson 
75fb90e4d1STodd Munson   PetscReal mu1_i; /*  used for model agreement in interpolation */
76fb90e4d1STodd Munson   PetscReal mu2_i; /*  used for model agreement in interpolation */
77fb90e4d1STodd Munson 
78fb90e4d1STodd Munson   PetscReal gamma1_i; /*  factor used for interpolation */
79fb90e4d1STodd Munson   PetscReal gamma2_i; /*  factor used for interpolation */
80fb90e4d1STodd Munson   PetscReal gamma3_i; /*  factor used for interpolation */
81fb90e4d1STodd Munson   PetscReal gamma4_i; /*  factor used for interpolation */
82fb90e4d1STodd Munson 
83fb90e4d1STodd Munson   PetscReal theta_i; /*  factor used for interpolation */
84fb90e4d1STodd Munson 
85fb90e4d1STodd Munson   PetscReal min_radius; /*  lower bound on initial radius value */
86fb90e4d1STodd Munson   PetscReal max_radius; /*  upper bound on trust region radius */
87fb90e4d1STodd Munson   PetscReal epsilon;    /*  tolerance used when computing actred/prered */
88fb90e4d1STodd Munson 
89fb90e4d1STodd Munson   PetscInt init_type;   /*  Trust-region initialization method */
90fb90e4d1STodd Munson   PetscInt update_type; /*  Trust-region update method */
91fb90e4d1STodd Munson } TAO_NTR;
92