double NOX::LineSearch::Utils::Slope::
computeSlope(const Abstract::Vector& dir, const Abstract::Group& grp) 
{
   if (grp.isGradient()) 
     return(dir.innerProduct(grp.getGradient()));

  // Allocate space for vecPtr if necessary
   if (Teuchos::is_null(vecPtr)) 
     vecPtr = dir.clone(ShapeCopy);

  // v = J * dir
  NOX::Abstract::Group::ReturnType status = grp.applyJacobian(dir,*vecPtr);
  
  if (status != NOX::Abstract::Group::Ok) 
  {
    utils.out() << "NOX::LineSearch::Utils::Slope::computeSlope -  Unable to apply Jacobian!" << std::endl;
    throw "NOX Error";
  }

  // Check that F exists
  if (!grp.isF()) 
  {
    utils.out() << "NOX::LineSearch::Utils::Slope::computeSlope - Invalid F" << std::endl;
    throw "NOX Error";
  }

  // Return <v, F> = F' * J * dir = <J'F, dir> = <g, dir>
  return(vecPtr->innerProduct(grp.getF()));
}
Esempio n. 2
0
bool NOX::Direction::SteepestDescent::compute(Abstract::Vector& dir, 
				 Abstract::Group& soln, 
				 const Solver::Generic& solver) 
{
  NOX::Abstract::Group::ReturnType status;

  // Compute F at current solution
  status = soln.computeF();
  if (status != NOX::Abstract::Group::Ok) 
    throwError("compute", "Unable to compute F");

  // Compute Jacobian at current solution
  status = soln.computeJacobian();
  if (status != NOX::Abstract::Group::Ok) 
    throwError("compute", "Unable to compute Jacobian");

  // Scale
  switch (scaleType) {

  case NOX::Direction::SteepestDescent::TwoNorm:
    
    meritFuncPtr->computeGradient(soln, dir);
    dir.scale(-1.0/dir.norm());
    break;
    
  case NOX::Direction::SteepestDescent::FunctionTwoNorm:
    
    meritFuncPtr->computeGradient(soln, dir);
    dir.scale(-1.0/soln.getNormF());
    break;
    
  case NOX::Direction::SteepestDescent::QuadMin:
    
    meritFuncPtr->computeQuadraticMinimizer(soln, dir);
      
    break;
    
  case NOX::Direction::SteepestDescent::None:
    
    meritFuncPtr->computeGradient(soln, dir);
    dir.scale( -1.0 );
    break;

  default:
    
    throwError("compute", "Invalid scaleType");
    
  }

  return true;
}
bool NOX::LineSearch::SafeguardedStep::compute(Abstract::Group& newGrp,
                           double& step,
                           const Abstract::Vector& dir,
                           const Solver::Generic& s)
{
  printOpeningRemarks();

  if (useCounter_) {
    counter_.incrementNumLineSearches();
    counter_.incrementNumNonTrivialLineSearches();
  }

  invLimits_->reciprocal(*userLimits_);
  *scaledUpdate_ = dir;
  scaledUpdate_->scale(*invLimits_);
  double infNorm = scaledUpdate_->norm(NOX::Abstract::Vector::MaxNorm);

  if (infNorm > 1.0)
    step = 1.0 / infNorm;
  else
    step = 1.0;

  if (print_.isPrintType(NOX::Utils::Details)) {
    print_.out () << "userLimits_:" << std::endl;
    userLimits_->print(print_.out());
    print_.out () << "scaledUpdate_:" << std::endl;
    scaledUpdate_->print(print_.out());
    print_.out () << "dir:" << std::endl;
    dir.print(print_.out());
    print_.out () << "Old Solution:" << std::endl;
    newGrp.getX().print(print_.out());
  }

  if (print_.isPrintType(NOX::Utils::InnerIteration)) {
    print_.out () << "    computed step = " << step << std::endl;
  }

  step = std::max(step,lowerStepBound_);
  step = std::min(step,upperStepBound_);

  if (print_.isPrintType(NOX::Utils::InnerIteration))
    print_.out () << "    bounded step = " << step << std::endl;

  newGrp.computeX(newGrp,dir,step);

  if (print_.isPrintType(NOX::Utils::Details)) {
    // reuse invLimits_, will be reset above
    invLimits_->update(step,dir,0.0);
    print_.out () << "Final Step Scaled Update:" << std::endl;
    invLimits_->print(print_.out());
    print_.out () << "New Solution:" << std::endl;
    newGrp.getX().print(print_.out());
  }

  if (useCounter_)
    counter_.setValues(*paramsPtr_);

  return true;
}
double NOX::LineSearch::Utils::Slope::
computeSlopeWithOutJac(const Abstract::Vector& dir, 
		       const Abstract::Group& grp) 
{
  // Allocate space for vecPtr and grpPtr if necessary
  if (Teuchos::is_null(vecPtr)) 
    vecPtr = dir.clone(ShapeCopy);
  if (Teuchos::is_null(grpPtr))
    grpPtr = grp.clone(ShapeCopy);

  // Check that F exists
  if (!grp.isF()) 
  {
    utils.out() << "NOX::LineSearch::Utils::Slope::computeSlope - Invalid F" << std::endl;
    throw "NOX Error";
  }

  // Compute the perturbation parameter
  double lambda = 1.0e-6;
  double denominator = dir.norm();

  // Don't divide by zero
  if (denominator == 0.0)
    denominator = 1.0;

  double eta = lambda * (lambda + grp.getX().norm() / denominator);

  // Don't divide by zero
  if (eta == 0.0)
    eta = 1.0e-6;

  // Perturb the solution vector
  vecPtr->update(eta, dir, 1.0, grp.getX(), 0.0);

  // Compute the new F --> F(x + eta * dir)
  grpPtr->setX(*vecPtr);  
  grpPtr->computeF();

  // Compute Js = (F(x + eta * dir) - F(x))/eta
  vecPtr->update(-1.0/eta, grp.getF(), 1.0/eta, grpPtr->getF(), 0.0);
  
  return(vecPtr->innerProduct(grp.getF()));
}
Esempio n. 5
0
void NOX::Direction::QuasiNewton::MemoryUnit::reset(const Abstract::Vector& newX, 
							 const Abstract::Vector& oldX, 
							 const Abstract::Vector& newG, 
							 const Abstract::Vector& oldG)
{
  if (Teuchos::is_null(sPtr))
  {
    sPtr = newX.clone(ShapeCopy);
    yPtr = newX.clone(ShapeCopy);
  }

  sPtr->update(1.0, newX, -1.0, oldX, 0.0);
  yPtr->update(1.0, newG, -1.0, oldG, 0.0);
  sdotyValue = sPtr->innerProduct(*yPtr);
  ydotyValue = yPtr->innerProduct(*yPtr);
  rhoValue = 1.0 / sdotyValue;
}
bool NonlinearCG::compute(Abstract::Vector& dir, Abstract::Group& soln,
                          const Solver::Generic& solver)
{
  Abstract::Group::ReturnType ok;

  // Initialize vector memory if haven't already
  if(Teuchos::is_null(oldDirPtr))
    oldDirPtr = soln.getX().clone(NOX::ShapeCopy);
  if(Teuchos::is_null(oldDescentDirPtr))
    oldDescentDirPtr = soln.getX().clone(NOX::ShapeCopy);
  // These are conditionally created
  if(Teuchos::is_null(diffVecPtr) && usePRbeta)
    diffVecPtr = soln.getX().clone(NOX::ShapeCopy);
  if(Teuchos::is_null(tmpVecPtr) && doPrecondition)
    tmpVecPtr = soln.getX().clone(NOX::ShapeCopy);

  // Get a reference to the old solution group (const)
  oldSolnPtr = &solver.getPreviousSolutionGroup();
  const Abstract::Group& oldSoln(*oldSolnPtr);

  niter = solver.getNumIterations();

  // Construct Residual and precondition (if desired) as first step in 
  // getting new search direction

  ok = soln.computeF();
  if (ok != Abstract::Group::Ok) 
  {
    if (utils->isPrintType(Utils::Warning))
      utils->out() << "NOX::Direction::NonlinearCG::compute - Unable to compute F." << std::endl;
    return false;
  }

  dir = soln.getF();  

  if(doPrecondition) 
  {
    if(!soln.isJacobian())
      ok = soln.computeJacobian();
    if (ok != Abstract::Group::Ok) 
    {
      if (utils->isPrintType(Utils::Warning))
        utils->out() << "NOX::Direction::NonlinearCG::compute - Unable to compute Jacobian." << std::endl;
      return false;
    }

    *tmpVecPtr = dir;

    ok = soln.applyRightPreconditioning(false, paramsPtr->sublist("Nonlinear CG").sublist("Linear Solver"), *tmpVecPtr, dir);
    if( ok != Abstract::Group::Ok ) 
    {
      if (utils->isPrintType(Utils::Warning))
        utils->out() << "NOX::Direction::NonlinearCG::compute - Unable to apply Right Preconditioner." << std::endl;
      return false;
    }
  }

  dir.scale(-1.0);

  // Orthogonalize using previous search direction

  beta = 0.0;

  if( niter!=0 )
  {  
    // Two choices (for now) for orthogonalizing descent direction with previous:
    if( usePRbeta )
    {
      // Polak-Ribiere beta
      *diffVecPtr = dir;
      diffVecPtr->update(-1.0, *oldDescentDirPtr, 1.0); 

      double denominator = oldDescentDirPtr->innerProduct(oldSoln.getF());

      beta = diffVecPtr->innerProduct(soln.getF()) / denominator;

      // Constrain beta >= 0
      if( beta < 0.0 ) 
      {
        if (utils->isPrintType(Utils::OuterIteration))
          utils->out() << "BETA < 0, (" << beta << ") --> Resetting to zero" << std::endl;
        beta = 0.0;
      }
    } 
    else
    {
      // Fletcher-Reeves beta
      double denominator = oldDescentDirPtr->innerProduct(oldSoln.getF());

      beta = dir.innerProduct(soln.getF()) / denominator;

    } 

    //  Allow for restart after specified number of nonlinear iterations
    if( (niter % restartFrequency) == 0 )
    {
      if( utils->isPrintType(Utils::OuterIteration) )
        utils->out() << "Resetting beta --> 0" << std::endl;

      beta = 0 ;  // Restart with Steepest Descent direction
    }
  } // niter != 0

  *oldDescentDirPtr = dir;

  dir.update(beta, *oldDirPtr, 1.0);

  *oldDirPtr = dir;

  return (ok == Abstract::Group::Ok);
}