Actual source code: nls.h

  1: /*$Id$*/

  3: // Context for a Newton line search method (unconstrained minimization)

  5: #ifndef __TAO_NLS_H
  7: #include "tao_solver.h"
  8: #include "src/matrix/lmvmmat.h"

 10: typedef struct {
 11:   TaoLMVMMat *M;

 13:   TaoVec *G;
 14:   TaoVec *D;
 15:   TaoVec *W;

 17:   TaoVec *Xold;
 18:   TaoVec *Gold;
 19:   TaoVec *Diag;

 21:   // Parameters when updating the perturbation added to the Hessian matrix
 22:   double sval;          // Starting perturbation value, default zero
 23:                         
 24:   double imin;          // Minimum perturbation added during initialization 
 25:   double imax;          // Maximum perturbation added during initialization
 26:   double imfac;         // Merit function factor during initialization

 28:   double pmin;          // Minimim perturbation value
 29:   double pmax;          // Maximum perturbation value
 30:   double pgfac;         // Perturbation growth factor
 31:   double psfac;         // Perturbation shrink factor
 32:   double pmgfac;        // Merit function growth factor
 33:   double pmsfac;        // Merit function shrink factor

 35:   // The perturbation to the Hessian matrix is initialized and updated
 36:   // according to the following scheme:
 37:   //
 38:   //   pert = sval;
 39:   //
 40:   //   do until convergence
 41:   //     shift Hessian by pert
 42:   //     solve Newton system
 43:   //
 44:   //     if (linear solver failed or did not compute a descent direction)
 45:   //       use steepest descent direction and increase perturbation
 46:   //
 47:   //       if (0 == pert)
 48:   //         initialize perturbation
 49:   //         pert = min(imax, max(imin, imfac * norm(G)))
 50:   //       else
 51:   //         increase perturbation
 52:   //         pert = min(pmax, max(pgfac * pert, pmgfac * norm(G)))
 53:   //       fi
 54:   //     else
 55:   //       use linear solver direction and decrease perturbation
 56:   //
 57:   //       pert = min(psfac * pert, pmsfac * norm(G))
 58:   //       if (pert < pmin)
 59:   //         pert = 0
 60:   //       fi
 61:   //     fi
 62:   //
 63:   //     perform line search
 64:   //     function and gradient evaluation
 65:   //     check convergence
 66:   //   od

 68:   // Parameters when updating the trust-region radius based on steplength
 69:   double nu1;                // used to compute trust-region radius
 70:   double nu2;                // used to compute trust-region radius
 71:   double nu3;                // used to compute trust-region radius
 72:   double nu4;                // used to compute trust-region radius

 74:   double omega1;        // factor used for trust-region update
 75:   double omega2;        // factor used for trust-region update
 76:   double omega3;        // factor used for trust-region update
 77:   double omega4;        // factor used for trust-region update
 78:   double omega5;        // factor used for trust-region update

 80:   // if   step < nu1                  (very bad step)
 81:   //   radius = omega1 * min(norm(d), radius)
 82:   // elif step < nu2                (bad step)
 83:   //   radius = omega2 * min(norm(d), radius)
 84:   // elif step < nu3                (okay step)
 85:   //   radius = omega3 * radius;
 86:   // elif step < nu4                (good step)
 87:   //   radius = max(omega4 * norm(d), radius)
 88:   // else                         (very good step)
 89:   //   radius = max(omega5 * norm(d), radius)
 90:   // fi
 91:  
 92:   // Parameters when updating the trust-region radius based on reduction
 93:   double eta1;                // used to compute trust-region radius
 94:   double eta2;                // used to compute trust-region radius
 95:   double eta3;                // used to compute trust-region radius
 96:   double eta4;                // used to compute trust-region radius

 98:   double alpha1;        // factor used for trust-region update
 99:   double alpha2;        // factor used for trust-region update
100:   double alpha3;        // factor used for trust-region update
101:   double alpha4;        // factor used for trust-region update
102:   double alpha5;        // factor used for trust-region update

104:   // kappa = ared / pred
105:   // if   kappa < eta1                 (very bad step)
106:   //   radius = alpha1 * min(norm(d), radius)
107:   // elif kappa < eta2                (bad step)
108:   //   radius = alpha2 * min(norm(d), radius)
109:   // elif kappa < eta3                (okay step)
110:   //   radius = alpha3 * radius;
111:   // elif kappa < eta4                (good step)
112:   //   radius = max(alpha4 * norm(d), radius)
113:   // else                         (very good step)
114:   //   radius = max(alpha5 * norm(d), radius)
115:   // fi
116:  
117:   // Parameters when updating the trust-region radius based on interpolation
118:   double mu1;                // used for model agreement in interpolation
119:   double mu2;                // used for model agreement in interpolation

121:   double gamma1;        // factor used for interpolation
122:   double gamma2;        // factor used for interpolation
123:   double gamma3;        // factor used for interpolation
124:   double gamma4;        // factor used for interpolation

126:   double theta;                // factor used for interpolation

128:   // kappa = ared / pred
129:   // if   kappa >= 1.0 - mu1        (very good step)
130:   //   choose tau in [gamma3, gamma4]
131:   //   radius = max(tau * norm(d), radius)
132:   // elif kappa >= 1.0 - mu2    (good step)
133:   //   choose tau in [gamma2, gamma3]
134:   //   if (tau >= 1.0)
135:   //     radius = max(tau * norm(d), radius)
136:   //   else
137:   //     radius = tau * min(norm(d), radius)
138:   //   fi
139:   // else                         (bad step)
140:   //   choose tau in [gamma1, 1.0]
141:   //   radius = tau * min(norm(d), radius)
142:   // fi
143:  
144:   // Parameters when initializing trust-region radius based on interpolation
145:   double mu1_i;                // used for model agreement in interpolation
146:   double mu2_i;                // used for model agreement in interpolation

148:   double gamma1_i;        // factor used for interpolation
149:   double gamma2_i;        // factor used for interpolation
150:   double gamma3_i;        // factor used for interpolation
151:   double gamma4_i;        // factor used for interpolation

153:   double theta_i;        // factor used for interpolation

155:   // Other parameters
156:   double min_radius;    // lower bound on initial radius value
157:   double max_radius;    // upper bound on trust region radius
158:   double epsilon;       // tolerance used when computing ared/pred

160:   int newt;                // Newton directions attempted
161:   int bfgs;                // BFGS directions attempted
162:   int sgrad;                // Scaled gradient directions attempted
163:   int grad;                // Gradient directions attempted

165:   int ksp_type;                // KSP method for the code
166:   int pc_type;                // Preconditioner for the code
167:   int bfgs_scale_type;        // Scaling matrix to used for the bfgs preconditioner
168:   int init_type;        // Trust-region initialization method
169:   int update_type;      // Trust-region update method
170: } TAO_NLS;

172: #endif