void NonGradient::optimize_vertex_positions(PatchData &pd, 
                                            MsqError &err)
{
  MSQ_FUNCTION_TIMER( "NonGradient::optimize_vertex_positions" );
  int numRow = getDimension();
  int numCol = numRow+1;  
  std::vector<double> height(numCol); 

  for(int col = 0; col < numCol; col++)
  {
    height[col] =  evaluate(&simplex[col*numRow], pd, err);//  eval patch stuff
  }
  if(mNonGradDebug > 0)
  {
    printSimplex( simplex, height );
  }

  // standardization
  TerminationCriterion* term_crit=get_inner_termination_criterion();
  int maxNumEval = getMaxNumEval();
  double threshold = getThreshold();
//  double ftol = getTolerance();
  int ilo = 0;  //height[ilo]<=...
  int inhi = 0; //...<=height[inhi]<=
  int ihi = 0;  //<=height[ihi] 
//  double rtol = 2.*ftol;
  double ysave;
  double ytry;
  bool afterEvaluation = false;
  std::vector<double> rowSum(numRow);
  getRowSum( numRow, numCol, simplex, rowSum);
  while( !( term_crit->terminate() ) )
  {

    if(mNonGradDebug > 0)
    {
      printSimplex( simplex, height );
    }
    //std::cout << "rtol " << rtol << " ftol " << ftol << " MesquiteIter " << term_crit->get_iteration_count() << " Done " << term_crit->terminate() << std::endl;

    if( afterEvaluation )   
    {
      // reflect highPt through opposite face
      // height[0] may vanish
/*
      if( !testRowSum( numRow, numCol, &simplex[0], &rowSum[0]) )
      {
        // Before uncommenting here and ... 
        //MSQ_SETERR(err)("Internal check sum test A failed", MsqError::INTERNAL_ERROR);
        //MSQ_ERRRTN(err);
      }
*/
      ytry=amotry(simplex,height,&rowSum[0],ihi,-1.0, pd, err);
/*
      if( !testRowSum( numRow, numCol, &simplex[0], &rowSum[0]) )
      {
         // ... here, determine a maxVal majorizing the previous as well as the current simplex.
         //MSQ_SETERR(err)("Internal check sum test B failed", MsqError::INTERNAL_ERROR);
         //MSQ_ERRRTN(err);
      }   
*/
  
/*
      if( height[0] == 0.)
      {
         MSQ_SETERR(err)("(B) Zero objective function value", MsqError::INTERNAL_ERROR);
         exit(-1);
      }
*/
      if (ytry <= height[ilo])   
      {
        ytry=amotry(simplex,height,&rowSum[0],ihi,-2.0,pd,err);
        if( mNonGradDebug >= 3 ) 
        {      
         std::cout << "Reflect and Expand from highPt " << ytry << std::endl;
        }      
        //MSQ_PRINT(3)("Reflect and Expand from highPt : %e\n",ytry);
      }
      else 
      {
        if (ytry >= height[inhi]) 
        {
          ysave=height[ihi]; // Contract along highPt
          ytry=amotry(simplex,height,&rowSum[0],ihi,0.5,pd,err);
          if (ytry >= ysave)
          { // contract all directions toward lowPt
            for (int col=0;col<numCol;col++)
            {
              if (col != ilo)
              {
                for (int row=0;row<numRow;row++)
                {
                  rowSum[row]=0.5*(simplex[row+col*numRow]+simplex[row+ilo*numRow]);
                  simplex[row+col*numRow]=rowSum[row];
                }
                height[col] = evaluate(&rowSum[0], pd, err); 
                if( mNonGradDebug >= 3 ) 
                {      
                  std::cout << "Contract all directions toward lowPt value( " << col << " ) = " << height[col] << " ilo = " << ilo << std::endl;
                }      
                //MSQ_PRINT(3)("Contract all directions toward lowPt value( %d ) = %e    ilo = %d\n", col, height[col], ilo);
              }
            }
          }
        }   
      } // ytri > h(ilo) 
    } // after evaluation
    ilo=1; // conditional operator or inline if 
    ihi = height[0] > height[1] ? (inhi=1,0) : (inhi=0,1);
    for (int col=0;col<numCol;col++)
    {
      if (height[col] <= height[ilo])
      {
        ilo=col;  // ilo := argmin height
      }
      if (height[col] > height[ihi])
      {
        inhi=ihi;
        ihi=col;
      } 
      else  // height[ihi] >= height[col]
        if (col != ihi && height[col] > height[inhi] ) inhi=col;
    }
//    rtol=2.0*fabs( height[ihi]-height[ilo] )/
//         ( fabs(height[ihi])+fabs(height[ilo])+threshold );
    afterEvaluation = true;
  } //  while not converged 

  // Always set to current best mesh.
  { 
    if( ilo != 0 )
    {
      double yTemp = height[0];
      height[0] = height[ilo]; // height dimension numCol
      height[ilo] = yTemp;
      for (int row=1;row<numRow;row++)
      { 
          yTemp = simplex[row];
          simplex[row] = simplex[row+ilo*numRow];
          simplex[row+ilo*numRow] = yTemp;
      }
    }
  }
  if( pd.num_free_vertices() > 1 )
  {
    MSQ_SETERR(err)("Only one free vertex per patch implemented", MsqError::NOT_IMPLEMENTED);
  }

  Vector3D newPoint( &simplex[0] ); 
  size_t vertexIndex = 0; // fix c.f. freeVertexIndex
  pd.set_vertex_coordinates( newPoint, vertexIndex, err ); 
  pd.snap_vertex_to_domain( vertexIndex, err );
  if( term_crit->terminate() )
  {
    if( mNonGradDebug >= 1 ) 
    {      
         std::cout << "Optimization Termination OptStatus: Max Iter Exceeded" << std::endl;
    }      
    //MSQ_PRINT(1)("Optimization Termination OptStatus: Max Iter Exceeded\n"); 
  }
}
NonGradient::NonGradient(ObjectiveFunction* of, MsqError &err) 
  : VertexMover(of),
    PatchSetUser(true),
    projectGradient(false),
    mDimension(0),
    mThreshold(0.0),
    mTolerance(0.0),
    mMaxNumEval(0),
    mNonGradDebug(0),
    mUseExactPenaltyFunction(true),
    mScaleDiameter(0.1)
{
  set_debugging_level(2);
  //set the default inner termination criterion
  TerminationCriterion* default_crit=get_inner_termination_criterion();
  if(default_crit==NULL){
    MSQ_SETERR(err)("QualityImprover did not create a default inner "
                    "termination criterion.", MsqError::INVALID_STATE);
    return;
  }
  else{
    default_crit->add_iteration_limit( 5 );
    MSQ_ERRRTN(err);
  }
}  
double 
NonGradient::evaluate( double *point,  PatchData &pd, MsqError &err )
{
  if( pd.num_free_vertices() > 1 )
  {
    MSQ_SETERR(err)("Only one free vertex per patch implemented", MsqError::NOT_IMPLEMENTED);
  }
  const size_t vertexIndex = 0; 
  Vector3D originalVec = pd.vertex_by_index(vertexIndex);
  Vector3D pointVec;
  for( int index = 0; index<3; index++)
  {
    pointVec[index] = point[index];
  }
  pd.set_vertex_coordinates( pointVec, vertexIndex, err ); 
  pd.snap_vertex_to_domain( vertexIndex, err );  

  OFEvaluator& obj_func = get_objective_function_evaluator();
  TerminationCriterion* term_crit=get_inner_termination_criterion();

  double value;
  bool feasible = obj_func.evaluate( pd, value, err ); // MSQ_ERRZERO(err);
  term_crit->accumulate_inner( pd, value, 0, err );  //MSQ_CHKERR(err);
  if (err.error_code() == err.BARRIER_VIOLATED)
    err.clear();   // barrier violation not an error here
  MSQ_ERRZERO(err);
  pd.set_vertex_coordinates( originalVec, vertexIndex, err ); 
  if( !feasible && mUseExactPenaltyFunction )  
  { // "value" undefined btw
    double ensureFiniteRtol= .25;
    value = DBL_MAX * ensureFiniteRtol;
    //std::cout << " Infeasible patch: " << value << std::endl; printPatch( pd, err );
  }
  return value;
}
Exemplo n.º 4
0
ConjugateGradient::ConjugateGradient(ObjectiveFunction* objective,
                                     MsqError &err) :
  VertexMover(objective),
  PatchSetUser(true),
  pMemento(NULL),
  conjGradDebug(0)
{
    //Michael:: default to global?
  set_debugging_level(0);
    //set the default inner termination criterion
  TerminationCriterion* default_crit=get_inner_termination_criterion();
  if(default_crit==NULL){
    MSQ_SETERR(err)("QualityImprover did not create a default inner "
                    "termination criterion.", MsqError::INVALID_STATE);
    return;
  }
  else{
    default_crit->add_iteration_limit( 5 );
    MSQ_ERRRTN(err);
  }
  
}  
NonGradient::NonGradient(ObjectiveFunction* of) 
  : VertexMover(of),
    PatchSetUser(true),
    projectGradient(false),
    mDimension(0),
    mThreshold(0.0),
    mTolerance(0.0),
    mMaxNumEval(0),
    mNonGradDebug(0),
    mUseExactPenaltyFunction(true),
    mScaleDiameter(0.1)
{
  set_debugging_level(2);
  //set the default inner termination criterion
  TerminationCriterion* default_crit=get_inner_termination_criterion();
  if(default_crit==NULL){
    return;
  }
  else{
    default_crit->add_iteration_limit( 5 );
  }
}  
Exemplo n.º 6
0
/*!Performs Conjugate gradient minimization on the PatchData, pd.*/
void ConjugateGradient::optimize_vertex_positions(PatchData &pd, 
                                                MsqError &err){
  // pd.reorder();

  MSQ_FUNCTION_TIMER( "ConjugateGradient::optimize_vertex_positions" );

  Timer c_timer;
  size_t num_vert=pd.num_free_vertices();
  if(num_vert<1){
     MSQ_DBGOUT(1) << "\nEmpty free vertex list in ConjugateGradient\n";
     return;
  }
/*  
    //zero out arrays
  int zero_loop=0;
  while(zero_loop<arraySize){
    fGrad[zero_loop].set(0,0,0);
    pGrad[zero_loop].set(0,0,0);
    fNewGrad[zero_loop].set(0,0,0);
    ++zero_loop;
  }
*/
  
  // get OF evaluator
  OFEvaluator& objFunc = get_objective_function_evaluator();
  
  size_t ind;
    //Michael cull list:  possibly set soft_fixed flags here
  
   //MsqFreeVertexIndexIterator free_iter(pd, err);  MSQ_ERRRTN(err);
  
      
   double f=0;
   //Michael, this isn't equivalent to CUBIT because we only want to check
   //the objective function value of the 'bad' elements
   //if invalid initial patch set an error.
   bool temp_bool = objFunc.update(pd, f, fGrad, err);
   assert(fGrad.size() == num_vert);
   if(MSQ_CHKERR(err))
      return;
  if( ! temp_bool){
    MSQ_SETERR(err)("Conjugate Gradient not able to get valid gradient "
                    "and function values on intial patch.", 
                    MsqError::INVALID_MESH);
    return;
  }
  double grad_norm=MSQ_MAX_CAP;
  
  if(conjGradDebug>0){
    MSQ_PRINT(2)("\nCG's DEGUB LEVEL = %i \n",conjGradDebug);
    grad_norm=Linf(arrptr(fGrad),fGrad.size());
    MSQ_PRINT(2)("\nCG's FIRST VALUE = %f,grad_norm = %f",f,grad_norm);
    MSQ_PRINT(2)("\n   TIME %f",c_timer.since_birth());
    grad_norm=MSQ_MAX_CAP;
  }

    //Initializing pGrad (search direction).  
  pGrad.resize(fGrad.size());
  for (ind = 0; ind < num_vert; ++ind)
    pGrad[ind]=(-fGrad[ind]);

  int j=0; // total nb of step size changes ... not used much
  int i=0; // iteration counter
  unsigned m=0; // 
  double alp=MSQ_MAX_CAP; // alp: scale factor of search direction
    //we know inner_criterion is false because it was checked in
    //loop_over_mesh before being sent here.
  TerminationCriterion* term_crit=get_inner_termination_criterion();
  
    //while ((i<maxIteration && alp>stepBound && grad_norm>normGradientBound)
    //     && !inner_criterion){
  while(!term_crit->terminate()){
    ++i;
      //std::cout<<"\Michael delete i = "<<i;
    int k=0;
    alp=get_step(pd,f,k,err);
    j+=k;
    if(conjGradDebug>2){
      MSQ_PRINT(2)("\n  Alp initial, alp = %20.18f",alp);
    }

     // if alp == 0, revert to steepest descent search direction
    if(alp==0){
      for (m = 0; m < num_vert; ++m) {
        pGrad[m]=(-fGrad[m]);
      }
      alp=get_step(pd,f,k,err);
      j+=k;
      if(conjGradDebug>1){
        MSQ_PRINT(2)("\n CG's search direction reset.");
        if(conjGradDebug>2)
          MSQ_PRINT(2)("\n  Alp was zero, alp = %20.18f",alp);
      }
      
    }
    if(alp!=0){
      pd.move_free_vertices_constrained( arrptr(pGrad), num_vert, alp, err );
      MSQ_ERRRTN(err);
      
      if (! objFunc.update(pd, f, fNewGrad, err)){
        MSQ_SETERR(err)("Error inside Conjugate Gradient, vertices moved "
                        "making function value invalid.", 
                        MsqError::INVALID_MESH);
        return;
      }
      assert(fNewGrad.size() == (unsigned)num_vert);
      
      if(conjGradDebug>0){
        grad_norm=Linf(arrptr(fNewGrad),num_vert);
        MSQ_PRINT(2)("\nCG's VALUE = %f,  iter. = %i,  grad_norm = %f,  alp = %f",f,i,grad_norm,alp);
        MSQ_PRINT(2)("\n   TIME %f",c_timer.since_birth());
      }
      double s11=0;
      double s12=0;
      double s22=0;
      //free_iter.reset();
      //while (free_iter.next()) {
      //  m=free_iter.value();
      for (m = 0; m < num_vert; ++m) {
        s11+=fGrad[m]%fGrad[m];
        s12+=fGrad[m]%fNewGrad[m];
        s22+=fNewGrad[m]%fNewGrad[m];
      }
      
        // Steepest Descent (takes 2-3 times as long as P-R)
        //double bet=0;
      
        // Fletcher-Reeves (takes twice as long as P-R)
        //double bet = s22/s11;

        // Polack-Ribiere  
      double bet;  
      if (!divide( s22-s12, s11, bet ))
        return; // gradient is zero    
      //free_iter.reset();
      //while (free_iter.next()) {
      //  m=free_iter.value();
      for (m = 0; m < num_vert; ++m) {
        pGrad[m]=(-fNewGrad[m]+(bet*pGrad[m]));
        fGrad[m]=fNewGrad[m];
      }
      if(conjGradDebug>2){
        MSQ_PRINT(2)(" \nSEARCH DIRECTION INFINITY NORM = %e",
                   Linf(arrptr(fNewGrad),num_vert));
      }
      
    }//end if on alp == 0

    term_crit->accumulate_patch( pd, err ); MSQ_ERRRTN(err);
    term_crit->accumulate_inner( pd, f, arrptr(fGrad), err );  MSQ_ERRRTN(err);
  }//end while
  if(conjGradDebug>0){
    MSQ_PRINT(2)("\nConjugate Gradient complete i=%i ",i);
    MSQ_PRINT(2)("\n-  FINAL value = %f, alp=%4.2e grad_norm=%4.2e",f,alp,grad_norm);
    MSQ_PRINT(2)("\n   FINAL TIME %f",c_timer.since_birth());
  }
}
Exemplo n.º 7
0
void QuasiNewton::optimize_vertex_positions( PatchData& pd, MsqError& err )
{
  TerminationCriterion& term = *get_inner_termination_criterion();
  OFEvaluator& func = get_objective_function_evaluator();
  
  const double sigma = 1e-4;
  const double beta0 = 0.25;
  const double beta1 = 0.80;
  const double tol1 = 1e-8;
  const double epsilon = 1e-10;

  double norm_r; //, norm_g;
  double alpha, beta;
  double obj, objn;

  size_t i;
  
    // Initialize stuff
  const size_t nn = pd.num_free_vertices();
  double a[QNVEC], b[QNVEC], r[QNVEC];
  for (i = 0; i < QNVEC; ++i)
    r[i] = 0;
  for (i = 0; i <= QNVEC; ++i) {
    v[i].clear();
    v[i].resize( nn, Vector3D(0.0) );
    w[i].clear();
    w[i].resize( nn, Vector3D(0.0) );
  }
  d.resize( nn );
  mHess.resize( nn );  //hMesh(mesh);

  bool valid = func.update( pd, obj, v[QNVEC], mHess, err ); MSQ_ERRRTN(err);
  if (!valid) {
    MSQ_SETERR(err)("Initial objective function is not valid", MsqError::INVALID_MESH);
    return;
  }

  while (!term.terminate()) {
    pd.recreate_vertices_memento( mMemento, err ); MSQ_ERRRTN(err);
    pd.get_free_vertex_coordinates( w[QNVEC] );

    x = v[QNVEC];
    for (i = QNVEC; i--; ) {
      a[i] = r[i] * inner( &(w[i][0]), arrptr(x), nn );
      plus_eq_scaled( arrptr(x), -a[i], &v[i][0], nn );
    }
     
    solve( arrptr(d), arrptr(x) );
  
    for (i = QNVEC; i--; ) {
      b[i] = r[i] * inner( &(v[i][0]), arrptr(d), nn );
      plus_eq_scaled( arrptr(d), a[i]-b[i], &(w[i][0]), nn );
    }
    
    alpha = -inner( &(v[QNVEC][0]), arrptr(d), nn );  /* direction is negated */
    if (alpha > 0.0) {
      MSQ_SETERR(err)("No descent.", MsqError::INVALID_MESH);
      return;
    }
   
    alpha *= sigma;
    beta = 1.0;
    
    pd.move_free_vertices_constrained( arrptr(d), nn, -beta, err ); MSQ_ERRRTN(err);
    valid = func.evaluate( pd, objn, v[QNVEC], err ); 
    if (err.error_code() == err.BARRIER_VIOLATED)             
      err.clear();  // barrier violated does not represent an actual error here
    MSQ_ERRRTN(err);
    if (!valid ||
        (obj - objn < -alpha*beta - epsilon &&
         length( &(v[QNVEC][0]), nn ) >= tol1)) {
      
      if (!valid)  // function not defined at trial point
        beta *= beta0;
      else  // unacceptable iterate
        beta *= beta1;
      
      for (;;) {
        if (beta < tol1) {
          pd.set_to_vertices_memento( mMemento, err ); MSQ_ERRRTN(err);
          MSQ_SETERR(err)("Newton step not good", MsqError::INTERNAL_ERROR);
          return;
        }
      
        pd.set_free_vertices_constrained( mMemento, arrptr(d), nn, -beta, err ); MSQ_ERRRTN(err);
        valid = func.evaluate( pd, objn, err );
        if (err.error_code() == err.BARRIER_VIOLATED)             
          err.clear();  // barrier violated does not represent an actual error here
        MSQ_ERRRTN(err);
        if (!valid) // function undefined at trial point
          beta *= beta0;
        else if (obj - objn < -alpha*beta - epsilon) // unacceptlable iterate
          beta *= beta1;
        else
          break;
      }
    }
    
    for (i = 0; i < QNVEC-1; ++i) {
      r[i] = r[i+1];
      w[i].swap( w[i+1] );
      v[i].swap( v[i+1] );
    }
    w[QNVEC-1].swap( w[0] );
    v[QNVEC-1].swap( v[0] );
    
    func.update( pd, obj, v[QNVEC], mHess, err ); MSQ_ERRRTN(err);
    norm_r = length_squared( &(v[QNVEC][0]), nn );
    //norm_g = sqrt(norm_r);

    // checks stopping criterion 
    term.accumulate_patch( pd, err ); MSQ_ERRRTN(err);
    term.accumulate_inner( pd, objn, &v[QNVEC][0], err ); MSQ_ERRRTN(err);
  }
}
Exemplo n.º 8
0
void TrustRegion::optimize_vertex_positions( PatchData& pd, MsqError& err )
{
  TerminationCriterion& term = *get_inner_termination_criterion();
  OFEvaluator& func = get_objective_function_evaluator();
  
  const double cg_tol = 1e-2;
  const double eta_1  = 0.01;
  const double eta_2  = 0.90;
  const double tr_incr = 10;
  const double tr_decr_def = 0.25;
  const double tr_decr_undef = 0.25;
  const double tr_num_tol = 1e-6;
  const int max_cg_iter = 10000;
  
  double radius = 1000;		/* delta*delta */

 const int nn = pd.num_free_vertices();
  wVect.resize(nn); Vector3D* w = arrptr(wVect);
  zVect.resize(nn); Vector3D* z = arrptr(zVect);
  dVect.resize(nn); Vector3D* d = arrptr(dVect);
  pVect.resize(nn); Vector3D* p = arrptr(pVect);
  rVect.resize(nn); Vector3D* r = arrptr(rVect);

  double norm_r, norm_g;
  double alpha, beta, kappa;
  double rz, rzm1;
  double dMp, norm_d, norm_dp1, norm_p;
  double obj, objn;

  int cg_iter;
  bool valid;

  mHess.initialize( pd, err );  //hMesh(mesh);
  valid = func.update( pd, obj, mGrad, mHess, err ); MSQ_ERRRTN(err);
  if (!valid) {
    MSQ_SETERR(err)("Initial objective function is not valid", MsqError::INVALID_MESH);
    return;
  }
  compute_preconditioner( err ); MSQ_ERRRTN(err);
  pd.recreate_vertices_memento( mMemento, err ); MSQ_ERRRTN(err);

  while (!term.terminate() && (radius > 1e-20)) {

    norm_r = length_squared(arrptr(mGrad), nn);
    norm_g = sqrt(norm_r);

    memset(d, 0, 3*sizeof(double)*nn);
    memcpy(r, arrptr(mGrad), nn*sizeof(Vector3D)); //memcpy(r, mesh->g, 3*sizeof(double)*nn);
    norm_g *= cg_tol;

    apply_preconditioner( z, r, err); MSQ_ERRRTN(err); //prec->apply(z, r, prec, mesh);
    negate(p, z, nn);
    rz = inner(r, z, nn);

    dMp    = 0;
    norm_p = rz;
    norm_d = 0;

    cg_iter = 0;
    while ((sqrt(norm_r) > norm_g) && 
#ifdef DO_STEEP_DESC
         (norm_d > tr_num_tol) && 
#endif
         (cg_iter < max_cg_iter)) 
    {
      ++cg_iter;

      memset(w, 0, 3*sizeof(double)*nn);
      //matmul(w, mHess, p); //matmul(w, mesh, p);
      mHess.product( w, p );

      kappa = inner(p, w, nn);
      if (kappa <= 0.0) {
        alpha = (sqrt(dMp*dMp+norm_p*(radius-norm_d))-dMp)/norm_p;
        plus_eq_scaled( d, alpha, p, nn );
	break;
      }

      alpha = rz / kappa;

      norm_dp1 = norm_d + 2.0*alpha*dMp + alpha*alpha*norm_p;
      if (norm_dp1 >= radius) {
        alpha = (sqrt(dMp*dMp+norm_p*(radius-norm_d))-dMp)/norm_p;
        plus_eq_scaled( d, alpha, p, nn );
	break;
      }

      plus_eq_scaled( d, alpha, p, nn );
      plus_eq_scaled( r, alpha, w, nn );
      norm_r = length_squared(r, nn);

      apply_preconditioner( z, r, err); MSQ_ERRRTN(err); //prec->apply(z, r, prec, mesh);

      rzm1 = rz;
      rz = inner(r, z, nn);
      beta = rz / rzm1;
      times_eq_minus( p, beta, z, nn );

      dMp = beta*(dMp + alpha*norm_p);
      norm_p = rz + beta*beta*norm_p;
      norm_d = norm_dp1;
    }

#ifdef DO_STEEP_DESC    
    if (norm_d <= tr_num_tol) {
      norm_g = length(arrptr(mGrad), nn);
      double ll = 1.0;
      if (norm_g < tr_num_tol)
        break;
      if (norm_g > radius)
        ll = radius / nurm_g;
      for (int i = 0; i < nn; ++i)
        d[i] = ll * mGrad[i];
    }
#endif

    alpha = inner( arrptr(mGrad), d, nn ); // inner(mesh->g, d, nn);

    memset(p, 0, 3*sizeof(double)*nn);
    //matmul(p, mHess, d); //matmul(p, mesh, d);
    mHess.product( p, d );
    beta = 0.5*inner(p, d, nn);
    kappa = alpha + beta;

    /* Put the new point into the locations */
    pd.move_free_vertices_constrained( d, nn, 1.0, err ); MSQ_ERRRTN(err);

    valid = func.evaluate( pd, objn, err ); MSQ_ERRRTN(err); 
    if (!valid) {
      /* Function not defined at trial point */
      radius *= tr_decr_undef;
      pd.set_to_vertices_memento( mMemento, err ); MSQ_ERRRTN(err); 
      continue;
    }
      

    if ((fabs(kappa) <= tr_num_tol) && (fabs(objn - obj) <= tr_num_tol)) {
      kappa = 1;
    }
    else {
      kappa = (objn - obj) / kappa;
    }
    
    if (kappa < eta_1) {
      /* Iterate is unacceptable */
      radius *= tr_decr_def;
      pd.set_to_vertices_memento( mMemento, err ); MSQ_ERRRTN(err); 
      continue;
    }

      /* Iterate is acceptable */
    if (kappa >= eta_2) {
      /* Iterate is a very good step, increase radius */
      radius *= tr_incr;
      if (radius > 1e20) {
	radius = 1e20;
      }
    }

    func.update( pd, obj, mGrad, mHess, err );
    compute_preconditioner( err ); MSQ_ERRRTN(err);
    pd.recreate_vertices_memento( mMemento, err ); MSQ_ERRRTN(err);

    // checks stopping criterion 
    term.accumulate_patch( pd, err ); MSQ_ERRRTN(err);
    term.accumulate_inner( pd, objn, arrptr(mGrad), err ); MSQ_ERRRTN(err);
  }
}    
void SteepestDescent::optimize_vertex_positions(PatchData &pd, 
                                                MsqError &err)
{
  MSQ_FUNCTION_TIMER( "SteepestDescent::optimize_vertex_positions" );

  const int SEARCH_MAX = 100;
  const double c1 = 1e-4;
  //std::vector<Vector3D> unprojected(pd.num_free_vertices()); 
  std::vector<Vector3D> gradient(pd.num_free_vertices()); 
  bool feasible=true;//bool for OF values
  double min_edge_len, max_edge_len;
  double step_size=0, original_value=0, new_value=0;
  double norm_squared=0;
  PatchDataVerticesMemento* pd_previous_coords;
  TerminationCriterion* term_crit=get_inner_termination_criterion();
  OFEvaluator& obj_func = get_objective_function_evaluator();
  
    // get vertex memento so we can restore vertex coordinates for bad steps.
  pd_previous_coords = pd.create_vertices_memento( err ); MSQ_ERRRTN(err);
    // use auto_ptr to automatically delete memento when we exit this function
  std::auto_ptr<PatchDataVerticesMemento> memento_deleter( pd_previous_coords );

    // Evaluate objective function.
    //
    // Always use 'update' version when beginning optimization so that
    // if doing block coordinate descent the OF code knows the set of
    // vertices we are modifying during the optimziation (the subset
    // of the mesh contained in the current patch.)  This has to be
    // done up-front because typically an OF will just store the portion
    // of the OF value (e.g. the numeric contribution to the sum for an
    // averaging OF) for the initial patch.
  feasible = obj_func.update( pd, original_value, gradient, err ); MSQ_ERRRTN(err);
    // calculate gradient dotted with itself
  norm_squared = length_squared( gradient );
  
    //set an error if initial patch is invalid.
  if(!feasible){
    MSQ_SETERR(err)("SteepestDescent passed invalid initial patch.",
                    MsqError::INVALID_ARG);
    return;
  }

    // use edge length as an initial guess for for step size
  pd.get_minmax_edge_length( min_edge_len, max_edge_len );
  //step_size = max_edge_len / std::sqrt(norm_squared);
  //if (!finite(step_size))  // zero-length gradient
  //  return;
//  if (norm_squared < DBL_EPSILON)
//    return;
  if (norm_squared >= DBL_EPSILON)
    step_size = max_edge_len / std::sqrt(norm_squared) * pd.num_free_vertices();

    // The steepest descent loop...
    // We loop until the user-specified termination criteria are met.
  while (!term_crit->terminate()) {
    MSQ_DBGOUT(3) << "Iteration " << term_crit->get_iteration_count() << std::endl;
    MSQ_DBGOUT(3) << "  o  original_value: " << original_value << std::endl;
    MSQ_DBGOUT(3) << "  o  grad norm suqared: " << norm_squared << std::endl;

      // Save current vertex coords so that they can be restored if
      // the step was bad.
    pd.recreate_vertices_memento( pd_previous_coords, err ); MSQ_ERRRTN(err);

      // Reduce step size until it satisfies Armijo condition
    int counter = 0;
    for (;;) {
      if (++counter > SEARCH_MAX || step_size < DBL_EPSILON) {
        MSQ_DBGOUT(3) << "    o  No valid step found.  Giving Up." << std::endl;
        return;
      }
      
      // Move vertices to new positions.
      // Note: step direction is -gradient so we pass +gradient and 
      //       -step_size to achieve the same thing.
      pd.move_free_vertices_constrained( arrptr(gradient), gradient.size(), -step_size, err ); MSQ_ERRRTN(err);
      // Evaluate objective function for new vertices.  We call the
      // 'evaluate' form here because we aren't sure yet if we want to
      // keep these vertices.  Until we call 'update', we have the option
      // of reverting a block coordinate decent objective function's state
      // to that of the initial vertex coordinates.  However, for block
      // coordinate decent to work correctly, we will need to call an
      // 'update' form if we decide to keep the new vertex coordinates.
      feasible = obj_func.evaluate( pd, new_value, err ); 
      if (err.error_code() == err.BARRIER_VIOLATED) 
        err.clear();  // barrier violated does not represent an actual error here
      MSQ_ERRRTN(err);
      MSQ_DBGOUT(3) << "    o  step_size: " << step_size << std::endl;
      MSQ_DBGOUT(3) << "    o  new_value: " << new_value << std::endl;

      if (!feasible) {
        // OF value is invalid, decrease step_size a lot
        step_size *= 0.2;
      }
      else if (new_value > original_value - c1 * step_size * norm_squared) {
        // Armijo condition not met.
        step_size *= 0.5;
      }
      else {
        // Armijo condition met, stop
        break;
      }
      
      // undo previous step : restore vertex coordinates
      pd.set_to_vertices_memento( pd_previous_coords, err );  MSQ_ERRRTN(err);
    }
   
      // Re-evaluate objective function to get gradient.
      // Calling the 'update' form here incorporates the new vertex 
      // positions into the 'accumulated' value if we are doing a 
      // block coordinate descent optimization.
    obj_func.update(pd, original_value, gradient, err ); MSQ_ERRRTN(err);
    if (projectGradient) {
      //if (cosineStep) {
      //  unprojected = gradient;
      //  pd.project_gradient( gradient, err ); MSQ_ERRRTN(err);
      //  double dot = inner_product( arrptr(gradient), arrptr(unprojected), gradient.size() );
      //  double lensqr1 = length_squared( gradient );
      //  double lensqr2 = length_squared( unprojected );
      //  double cossqr = dot * dot / lensqr1 / lensqr2;
      //  step_size *= sqrt(cossqr);
      //}
      //else {
        pd.project_gradient( gradient, err ); MSQ_ERRRTN(err);
      //}      
    }
      
      // Update terination criterion for next iteration.
      // This is necessary for efficiency.  Some values can be adjusted
      // for each iteration so we don't need to re-caculate the value
      // over the entire mesh.
    term_crit->accumulate_patch( pd, err );  MSQ_ERRRTN(err);
    term_crit->accumulate_inner( pd, original_value, arrptr(gradient), err ); MSQ_ERRRTN(err); 
      
      // Calculate initial step size for next iteration using step size 
      // from this iteration
    step_size *= norm_squared;
    norm_squared = length_squared( gradient );
//    if (norm_squared < DBL_EPSILON)
//      break;
  if (norm_squared >= DBL_EPSILON)
    step_size /= norm_squared;
  }
}