Actual source code: ntrimpl.h

  1: /*
  2:   Context for a Newton trust region method (unconstrained minimization)
  3: */

  5: #pragma once
  6: #include <petsc/private/taoimpl.h>

  8: typedef struct {
  9:   Mat M;
 10:   PC  bfgs_pre;

 12:   Vec D;
 13:   Vec W;

 15:   PetscReal radius;

 17:   /* Parameters when updating the trust-region radius based on reduction

 19:      kappa = ared / pred
 20:      if   kappa < eta1          (very bad step)
 21:        radius = alpha1 * min(norm(d), radius)
 22:      elif kappa < eta2          (bad step)
 23:        radius = alpha2 * min(norm(d), radius)
 24:      elif kappa < eta3          (okay step)
 25:        radius = alpha3 * radius;
 26:      elif kappa < eta4          (good step)
 27:        radius = max(alpha4 * norm(d), radius)
 28:      else                       (very good step)
 29:        radius = max(alpha5 * norm(d), radius)
 30:      fi
 31:   */

 33:   PetscReal eta1; /*  used to compute trust-region radius */
 34:   PetscReal eta2; /*  used to compute trust-region radius */
 35:   PetscReal eta3; /*  used to compute trust-region radius */
 36:   PetscReal eta4; /*  used to compute trust-region radius */

 38:   PetscReal alpha1; /*  factor used for trust-region update */
 39:   PetscReal alpha2; /*  factor used for trust-region update */
 40:   PetscReal alpha3; /*  factor used for trust-region update */
 41:   PetscReal alpha4; /*  factor used for trust-region update */
 42:   PetscReal alpha5; /*  factor used for trust-region update */

 44:   /* Parameters when updating the trust-region radius based on interpolation

 46:      kappa = ared / pred
 47:      if   kappa >= 1.0 - mu1    (very good step)
 48:        choose tau in [gamma3, gamma4]
 49:        radius = max(tau * norm(d), radius)
 50:      elif kappa >= 1.0 - mu2    (good step)
 51:        choose tau in [gamma2, gamma3]
 52:        if (tau >= 1.0)
 53:          radius = max(tau * norm(d), radius)
 54:        else
 55:          radius = tau * min(norm(d), radius)
 56:        fi
 57:      else                       (bad step)
 58:        choose tau in [gamma1, 1.0]
 59:        radius = tau * min(norm(d), radius)
 60:      fi
 61:   */

 63:   PetscReal mu1; /*  used for model agreement in radius update */
 64:   PetscReal mu2; /*  used for model agreement in radius update */

 66:   PetscReal gamma1; /*  factor used for radius update */
 67:   PetscReal gamma2; /*  factor used for radius update */
 68:   PetscReal gamma3; /*  factor used for radius update */
 69:   PetscReal gamma4; /*  factor used for radius update */

 71:   PetscReal theta; /*  factor used for radius update */

 73:   /* Parameters when initializing trust-region radius based on interpolation */

 75:   PetscReal mu1_i; /*  used for model agreement in interpolation */
 76:   PetscReal mu2_i; /*  used for model agreement in interpolation */

 78:   PetscReal gamma1_i; /*  factor used for interpolation */
 79:   PetscReal gamma2_i; /*  factor used for interpolation */
 80:   PetscReal gamma3_i; /*  factor used for interpolation */
 81:   PetscReal gamma4_i; /*  factor used for interpolation */

 83:   PetscReal theta_i; /*  factor used for interpolation */

 85:   PetscReal min_radius; /*  lower bound on initial radius value */
 86:   PetscReal max_radius; /*  upper bound on trust region radius */
 87:   PetscReal epsilon;    /*  tolerance used when computing actred/prered */

 89:   PetscInt init_type;   /*  Trust-region initialization method */
 90:   PetscInt update_type; /*  Trust-region update method */
 91: } TAO_NTR;