void NonGradient::printPatch(const PatchData &pd, MsqError &err)
{
  if( mNonGradDebug==0 )
  {
    return;
  }
  const size_t numNode = pd.num_nodes();   //27,  27 what?
  MSQ_PRINT(3)("Number of Vertices: %d\n",(int)pd.num_nodes());
  const size_t numVert = pd.num_free_vertices(); // 1
  MSQ_PRINT(3)("Num Free = %d\n",(int)pd.num_free_vertices());
  const size_t numSlaveVert = pd.num_slave_vertices(); //0
  const size_t numCoin = pd.num_corners(); // 64
  const MsqVertex* coord = pd.get_vertex_array(err);
  MSQ_PRINT(3)("Number of Vertices: %d\n",(int)pd.num_nodes());

  std::cout << "Patch " << numNode << "  " << numVert << "  " << numSlaveVert << "  " << numCoin << std::endl;
  MSQ_PRINT(3)("");
  std::cout << "Coordinate ";
  std::cout << "         " << std::endl;
  for( size_t index = 0; index < numVert; index++ )
  {
    std::cout << coord[index][0] << "  " << coord[index][1] << "  " << coord[index][2] << std::endl;
  }
  //const size_t numElt = pd.num_elements();
  if( mNonGradDebug >= 3 ) 
  {      
         std::cout << "Number of Elements: " << pd.num_elements() << std::endl;
  }      
  MSQ_PRINT(3)("Number of Elements: %d\n",(int)pd.num_elements());
}
bool
NonGradient::testRowSum( int numRow, int numCol, double* matrix, double* oldRowSum)
{
  bool same = true;
  std::vector<double> rowSum(numRow,0.);
  double maxVal = 0.;
  for (int col=0;col<numCol;col++)    
  {
    for (int row=0;row<numRow;row++)
    {
      rowSum[row] += matrix[row+col*numRow];
      if( fabs( matrix[row+col*numRow] ) > maxVal )
      {
         maxVal = fabs( matrix[row+col*numRow] );
      }
    }
  }
  double machEps = 1.e-14 * static_cast<double>(numRow); // better to use system parameters
  double upperBound = machEps * maxVal + 1.e-15;
  for (int row=0;row<numRow;row++)
  {
    if( fabs( rowSum[row] -  oldRowSum[row]) >  upperBound ) 
    { 
         same = false;
         if( mNonGradDebug >= 2 )
         {
           std::cout << " Warning: NonGradient Row Sum " << row << " Test: value " << rowSum[row] << " Discrepancy " <<  rowSum[row] -  oldRowSum[row] << " maxVal " << maxVal << " numRow " << numRow << " numCol " << numCol <<std::endl;
         }
         MSQ_PRINT(2)("NonGradient Row Sum [%d] Test failure: value %22.15e  difference %22.15e \n", row, rowSum[row], rowSum[row] -  oldRowSum[row]); 
    }
  }
  return(same);
}
void NonGradient::initialize_mesh_iteration(PatchData &pd, MsqError &err)
{
  int elementDimension = getPatchDimension( pd, err );  // to do: react to error
  int dimension = elementDimension * pd.num_free_vertices();
  //printPatch( pd, err );
  setDimension(dimension);
  int maxNumEval = 100*dimension;  // 1. Make this a user parameter
  setMaxNumEval(maxNumEval);
  double threshold = 1.e-10; // avoid division by zero
  setThreshold(threshold);
  double minEdgeLen = 0.0;
  double maxEdgeLen = 0.0;
//  double ftol = 0.;
  if( dimension > 0 )
  {
    pd.get_minmax_edge_length( minEdgeLen, maxEdgeLen );
    //ftol = minEdgeLen * 1.e-4; // Turn off Amoeba convergence criterion
    if( mNonGradDebug >= 1 ) 
    {      
         std::cout << "minimum edge length " << minEdgeLen << " maximum edge length " << maxEdgeLen << std::endl;
    }      
    MSQ_PRINT(3)("minimum edge length %e    maximum edge length %e\n", minEdgeLen,  maxEdgeLen);
  }
//  setTolerance(ftol);
  int numRow = dimension;
  int numCol = numRow+1;  
  if( numRow*numCol <= simplex.max_size() )
  { 
    simplex.assign(numRow*numCol, 0.);  // guard against previous simplex value
    double scale = minEdgeLen * mScaleDiameter;; 
    const MsqVertex* coord = pd.get_vertex_array(err);
    if( pd.num_free_vertices() > 1 )
    {
      MSQ_SETERR(err)("Only one free vertex per patch implemented", MsqError::NOT_IMPLEMENTED);
    }
    size_t index = 0;
    for( int col = 0; col < numCol; col++ )
    {
      for (int row=0;row<numRow;row++)
      {
        simplex[ row + col*numRow ] = coord[index][row];
        if( row == col-1 )
        {
          simplex[ row + col*numRow ] += scale/ static_cast<double>(numCol);
        }
      }
    }
  }
  else
  {
    MSQ_SETERR(err)("Use patch with fewer free vertices", MsqError::OUT_OF_MEMORY);
    if( mNonGradDebug >= 1 ) 
    {      
      std::cout << "ERROR: Too many free vertices in patch" << std::endl;
    }      
    //MSQ_PRINT(1)("ERROR: Too many free vertices in patch\n");
  }
}
// Extrapolate by a factor of fac through the face of the simplex
// opposite the high point to a new point.  If the new point is 
// better, swap it with the high point.
double
NonGradient::amotry( std::vector<double>& simplex, 
                 std::vector<double>& height,
                 double psum[], int ihi, double fac, PatchData &pd, MsqError &err)
{
  int numRow = getDimension();
  int numCol = numRow + 1;
  std::vector<double> ptry(numRow); // does this make sense?
  double fac1=(1.0-fac)/static_cast<double>(numRow);
  double fac2=fac1-fac;
  for (int row=0;row<numRow;row++)
  {
    ptry[row]=psum[row]*fac1-simplex[row+ihi*numRow]*fac2;
  }
  if( mNonGradDebug >= 3 ) 
  {      
    std::cout << "Try ";
  }      
  MSQ_PRINT(3)("Try");
  double ytry = evaluate(&ptry[0], pd, err); // value at trial point
  if( mNonGradDebug >= 3 ) 
  {      
    std::cout << ytry << std::endl;
  }      
  MSQ_PRINT(3)("yTry");
  if (ytry < height[ihi]) // better than highest (worst)
  {
    height[ihi]=ytry;     // swap ihi and ytry
    for (int row=0;row<numRow;row++)
    {
      psum[row] += (ptry[row]-simplex[row+ihi*numRow]);
      simplex[row+ihi*numRow]=ptry[row];
    }
  }
  return ytry;
}
Beispiel #5
0
double ConjugateGradient::get_step(PatchData &pd,double f0,int &j,
                                   MsqError &err)
{
  // get OF evaluator
  OFEvaluator& objFunc = get_objective_function_evaluator();

  size_t num_vertices=pd.num_free_vertices();
    //initial guess for alp
  double alp=1.0;
  int jmax=100;
  double rho=0.5;
    //feasible=false implies the mesh is not in the feasible region
  bool feasible=false;
  int found=0;
    //f and fnew hold the objective function value
  double f=0;
  double fnew=0;
    //Counter to avoid infinitly scaling alp
  j=0;
  //save memento
  pd.recreate_vertices_memento(pMemento, err);
    //if we must check feasiblility
    //while step takes mesh into infeasible region and ...
  while (j<jmax && !feasible && alp>MSQ_MIN) {
    ++j;
    pd.set_free_vertices_constrained(pMemento,arrptr(pGrad),num_vertices,alp,err);
    feasible=objFunc.evaluate(pd,f,err); MSQ_ERRZERO(err);
      //if not feasible, try a smaller alp (take smaller step)
    if(!feasible){
      alp*=rho;
    }
  }//end while ...
  
    //if above while ended due to j>=jmax, no valid step was found.
  if(j>=jmax){
    MSQ_PRINT(2)("\nFeasible Point Not Found");
    return 0.0;
  }
    //Message::print_info("\nOriginal f %f, first new f = %f, alp = %f",f0,f,alp);
    //if new f is larger than original, our step was too large
  if(f>=f0){
    j=0;
    while (j<jmax && found == 0){
      ++j;
      alp *= rho;
      pd.set_free_vertices_constrained(pMemento,arrptr(pGrad),num_vertices,alp,err);
        //Get new obj value
        //if patch is now invalid, then the feasible region is  convex or
        //we have an error.  For now, we assume an error.
      if(! objFunc.evaluate(pd,f,err) ){
        MSQ_SETERR(err)("Non-convex feasiblility region found.",MsqError::INVALID_MESH);
      }
      pd.set_to_vertices_memento(pMemento,err);MSQ_ERRZERO(err);
        //if our step has now improved the objective function value
      if(f<f0){
        found=1;
      }
    }//   end while j less than jmax
      //Message::print_info("\nj = %d found = %d f = %20.18f f0 = %20.18f\n",j,found,f,f0);
      //if above ended because of j>=jmax, take no step
    if(found==0){
        //Message::print_info("alp = %10.8f, but returning zero\n",alp);
      alp=0.0; 
      return alp;
    }

    j=0;
      //while shrinking the step improves the objFunc value further,
      //scale alp down.  Return alp, when scaling once more would
      //no longer improve the objFunc value.  
    while(j<jmax){
      ++j;
      alp*=rho;
      //step alp in search direction from original positions
      pd.set_free_vertices_constrained(pMemento,arrptr(pGrad),num_vertices,alp,err);MSQ_ERRZERO(err);

        //get new objective function value
      if (! objFunc.evaluate(pd,fnew,err))
        MSQ_SETERR(err)("Non-convex feasiblility region found while "
                        "computing new f.",MsqError::INVALID_MESH);
      if(fnew<f){
        f=fnew;
      }
      else{
	//Reset the vertices to original position
	pd.set_to_vertices_memento(pMemento,err);MSQ_ERRZERO(err);
	alp/=rho;
	return alp;
      }
    }
    //Reset the vertices to original position and return alp
    pd.set_to_vertices_memento(pMemento,err);MSQ_ERRZERO(err);
    return alp;
  }
    //else our new f was already smaller than our original
  else{
    j=0;
      //check to see how large of step we can take
    while (j<jmax && found == 0) {
      ++j;
        //scale alp up (rho must be less than 1)
      alp /= rho;
      //step alp in search direction from original positions
      pd.set_free_vertices_constrained(pMemento,arrptr(pGrad),num_vertices,alp,err);MSQ_ERRZERO(err);

      feasible = objFunc.evaluate(pd,fnew, err);MSQ_ERRZERO(err);
      if ( ! feasible ){
         alp *= rho;
	 
	 //Reset the vertices to original position and return alp
	 pd.set_to_vertices_memento(pMemento,err);MSQ_ERRZERO(err);
         return alp;
      }
      if (fnew<f) { 
          f = fnew; 
       }
       else {
         found=1;
         alp *= rho;
       }
    }

    //Reset the vertices to original position and return alp
    pd.set_to_vertices_memento(pMemento,err);MSQ_ERRZERO(err);
    return alp;
  }
}
Beispiel #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());
  }
}