xref: /petsc/src/tao/unconstrained/impls/nls/nlsimpl.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
11daac38eSTodd Munson /*
21daac38eSTodd Munson Context for a Newton line search method (unconstrained minimization)
31daac38eSTodd Munson */
41daac38eSTodd Munson 
5*a4963045SJacob Faibussowitsch #pragma once
61daac38eSTodd Munson #include <petsc/private/taoimpl.h>
71daac38eSTodd Munson 
81daac38eSTodd Munson typedef struct {
91daac38eSTodd Munson   Mat M;
100c51296cSAlp Dener   PC  bfgs_pre;
111daac38eSTodd Munson 
121daac38eSTodd Munson   Vec D;
131daac38eSTodd Munson   Vec W;
141daac38eSTodd Munson 
151daac38eSTodd Munson   Vec Xold;
161daac38eSTodd Munson   Vec Gold;
171daac38eSTodd Munson 
181daac38eSTodd Munson   /* Parameters when updating the perturbation added to the Hessian matrix
191daac38eSTodd Munson      according to the following scheme:
201daac38eSTodd Munson 
211daac38eSTodd Munson      pert = sval;
221daac38eSTodd Munson 
231daac38eSTodd Munson      do until convergence
241daac38eSTodd Munson        shift Hessian by pert
251daac38eSTodd Munson        solve Newton system
261daac38eSTodd Munson 
271daac38eSTodd Munson        if (linear solver failed or did not compute a descent direction)
281daac38eSTodd Munson          use steepest descent direction and increase perturbation
291daac38eSTodd Munson 
301daac38eSTodd Munson          if (0 == pert)
311daac38eSTodd Munson            initialize perturbation
321daac38eSTodd Munson            pert = min(imax, max(imin, imfac * norm(G)))
331daac38eSTodd Munson          else
341daac38eSTodd Munson            increase perturbation
351daac38eSTodd Munson            pert = min(pmax, max(pgfac * pert, pmgfac * norm(G)))
361daac38eSTodd Munson          fi
371daac38eSTodd Munson        else
381daac38eSTodd Munson          use linear solver direction and decrease perturbation
391daac38eSTodd Munson 
401daac38eSTodd Munson          pert = min(psfac * pert, pmsfac * norm(G))
411daac38eSTodd Munson          if (pert < pmin)
421daac38eSTodd Munson            pert = 0
431daac38eSTodd Munson          fi
441daac38eSTodd Munson        fi
451daac38eSTodd Munson 
461daac38eSTodd Munson        perform line search
471daac38eSTodd Munson        function and gradient evaluation
481daac38eSTodd Munson        check convergence
491daac38eSTodd Munson      od
501daac38eSTodd Munson   */
511daac38eSTodd Munson   PetscReal sval; /*  Starting perturbation value, default zero */
521daac38eSTodd Munson 
531daac38eSTodd Munson   PetscReal imin;  /*  Minimum perturbation added during initialization  */
541daac38eSTodd Munson   PetscReal imax;  /*  Maximum perturbation added during initialization */
551daac38eSTodd Munson   PetscReal imfac; /*  Merit function factor during initialization */
561daac38eSTodd Munson 
5735cb6cd3SPierre Jolivet   PetscReal pmin;   /*  Minimum perturbation value */
581daac38eSTodd Munson   PetscReal pmax;   /*  Maximum perturbation value */
591daac38eSTodd Munson   PetscReal pgfac;  /*  Perturbation growth factor */
601daac38eSTodd Munson   PetscReal psfac;  /*  Perturbation shrink factor */
611daac38eSTodd Munson   PetscReal pmgfac; /*  Merit function growth factor */
621daac38eSTodd Munson   PetscReal pmsfac; /*  Merit function shrink factor */
631daac38eSTodd Munson 
641daac38eSTodd Munson   /* Parameters when updating the trust-region radius based on steplength
651daac38eSTodd Munson      if   step < nu1            (very bad step)
661daac38eSTodd Munson        radius = omega1 * min(norm(d), radius)
671daac38eSTodd Munson      elif step < nu2            (bad step)
681daac38eSTodd Munson        radius = omega2 * min(norm(d), radius)
691daac38eSTodd Munson      elif step < nu3            (okay step)
701daac38eSTodd Munson        radius = omega3 * radius;
711daac38eSTodd Munson      elif step < nu4            (good step)
721daac38eSTodd Munson        radius = max(omega4 * norm(d), radius)
731daac38eSTodd Munson      else                       (very good step)
741daac38eSTodd Munson        radius = max(omega5 * norm(d), radius)
751daac38eSTodd Munson      fi
761daac38eSTodd Munson   */
771daac38eSTodd Munson   PetscReal nu1; /*  used to compute trust-region radius */
781daac38eSTodd Munson   PetscReal nu2; /*  used to compute trust-region radius */
791daac38eSTodd Munson   PetscReal nu3; /*  used to compute trust-region radius */
801daac38eSTodd Munson   PetscReal nu4; /*  used to compute trust-region radius */
811daac38eSTodd Munson 
821daac38eSTodd Munson   PetscReal omega1; /*  factor used for trust-region update */
831daac38eSTodd Munson   PetscReal omega2; /*  factor used for trust-region update */
841daac38eSTodd Munson   PetscReal omega3; /*  factor used for trust-region update */
851daac38eSTodd Munson   PetscReal omega4; /*  factor used for trust-region update */
861daac38eSTodd Munson   PetscReal omega5; /*  factor used for trust-region update */
871daac38eSTodd Munson 
881daac38eSTodd Munson   /* Parameters when updating the trust-region radius based on reduction
891daac38eSTodd Munson 
901daac38eSTodd Munson      kappa = ared / pred
911daac38eSTodd Munson      if   kappa < eta1          (very bad step)
921daac38eSTodd Munson        radius = alpha1 * min(norm(d), radius)
931daac38eSTodd Munson      elif kappa < eta2          (bad step)
941daac38eSTodd Munson        radius = alpha2 * min(norm(d), radius)
951daac38eSTodd Munson      elif kappa < eta3          (okay step)
961daac38eSTodd Munson        radius = alpha3 * radius;
971daac38eSTodd Munson      elif kappa < eta4          (good step)
981daac38eSTodd Munson        radius = max(alpha4 * norm(d), radius)
991daac38eSTodd Munson      else                       (very good step)
1001daac38eSTodd Munson        radius = max(alpha5 * norm(d), radius)
1011daac38eSTodd Munson      fi
1021daac38eSTodd Munson   */
1031daac38eSTodd Munson   PetscReal eta1; /*  used to compute trust-region radius */
1041daac38eSTodd Munson   PetscReal eta2; /*  used to compute trust-region radius */
1051daac38eSTodd Munson   PetscReal eta3; /*  used to compute trust-region radius */
1061daac38eSTodd Munson   PetscReal eta4; /*  used to compute trust-region radius */
1071daac38eSTodd Munson 
1081daac38eSTodd Munson   PetscReal alpha1; /*  factor used for trust-region update */
1091daac38eSTodd Munson   PetscReal alpha2; /*  factor used for trust-region update */
1101daac38eSTodd Munson   PetscReal alpha3; /*  factor used for trust-region update */
1111daac38eSTodd Munson   PetscReal alpha4; /*  factor used for trust-region update */
1121daac38eSTodd Munson   PetscReal alpha5; /*  factor used for trust-region update */
1131daac38eSTodd Munson 
1141daac38eSTodd Munson   /* Parameters when updating the trust-region radius based on interpolation
1151daac38eSTodd Munson 
1161daac38eSTodd Munson      kappa = ared / pred
1171daac38eSTodd Munson      if   kappa >= 1.0 - mu1    (very good step)
1181daac38eSTodd Munson        choose tau in [gamma3, gamma4]
1191daac38eSTodd Munson        radius = max(tau * norm(d), radius)
1201daac38eSTodd Munson      elif kappa >= 1.0 - mu2    (good step)
1211daac38eSTodd Munson        choose tau in [gamma2, gamma3]
1221daac38eSTodd Munson        if (tau >= 1.0)
1231daac38eSTodd Munson          radius = max(tau * norm(d), radius)
1241daac38eSTodd Munson        else
1251daac38eSTodd Munson          radius = tau * min(norm(d), radius)
1261daac38eSTodd Munson        fi
1271daac38eSTodd Munson      else                       (bad step)
1281daac38eSTodd Munson        choose tau in [gamma1, 1.0]
1291daac38eSTodd Munson        radius = tau * min(norm(d), radius)
1301daac38eSTodd Munson      fi
1311daac38eSTodd Munson   */
1321daac38eSTodd Munson   PetscReal mu1; /*  used for model agreement in interpolation */
1331daac38eSTodd Munson   PetscReal mu2; /*  used for model agreement in interpolation */
1341daac38eSTodd Munson 
1351daac38eSTodd Munson   PetscReal gamma1; /*  factor used for interpolation */
1361daac38eSTodd Munson   PetscReal gamma2; /*  factor used for interpolation */
1371daac38eSTodd Munson   PetscReal gamma3; /*  factor used for interpolation */
1381daac38eSTodd Munson   PetscReal gamma4; /*  factor used for interpolation */
1391daac38eSTodd Munson 
1401daac38eSTodd Munson   PetscReal theta; /*  factor used for interpolation */
1411daac38eSTodd Munson 
1421daac38eSTodd Munson   /*  Parameters when initializing trust-region radius based on interpolation */
1431daac38eSTodd Munson   PetscReal mu1_i; /*  used for model agreement in interpolation */
1441daac38eSTodd Munson   PetscReal mu2_i; /*  used for model agreement in interpolation */
1451daac38eSTodd Munson 
1461daac38eSTodd Munson   PetscReal gamma1_i; /*  factor used for interpolation */
1471daac38eSTodd Munson   PetscReal gamma2_i; /*  factor used for interpolation */
1481daac38eSTodd Munson   PetscReal gamma3_i; /*  factor used for interpolation */
1491daac38eSTodd Munson   PetscReal gamma4_i; /*  factor used for interpolation */
1501daac38eSTodd Munson 
1511daac38eSTodd Munson   PetscReal theta_i; /*  factor used for interpolation */
1521daac38eSTodd Munson 
1531daac38eSTodd Munson   /*  Other parameters */
1541daac38eSTodd Munson   PetscReal min_radius; /*  lower bound on initial radius value */
1551daac38eSTodd Munson   PetscReal max_radius; /*  upper bound on trust region radius */
1561daac38eSTodd Munson   PetscReal epsilon;    /*  tolerance used when computing ared/pred */
1571daac38eSTodd Munson 
1581daac38eSTodd Munson   PetscInt newt; /*  Newton directions attempted */
1591daac38eSTodd Munson   PetscInt bfgs; /*  BFGS directions attempted */
1601daac38eSTodd Munson   PetscInt grad; /*  Gradient directions attempted */
1611daac38eSTodd Munson 
1621daac38eSTodd Munson   PetscInt init_type;   /*  Trust-region initialization method */
1631daac38eSTodd Munson   PetscInt update_type; /*  Trust-region update method */
1641daac38eSTodd Munson 
1651daac38eSTodd Munson   PetscInt ksp_atol;
1661daac38eSTodd Munson   PetscInt ksp_rtol;
1671daac38eSTodd Munson   PetscInt ksp_ctol;
1681daac38eSTodd Munson   PetscInt ksp_negc;
1691daac38eSTodd Munson   PetscInt ksp_dtol;
1701daac38eSTodd Munson   PetscInt ksp_iter;
1711daac38eSTodd Munson   PetscInt ksp_othr;
1721daac38eSTodd Munson } TAO_NLS;
173