bool
NOX::Solver::TensorBased::implementGlobalStrategy(NOX::Abstract::Group& newGrp,
                      double& in_stepSize,
                      const NOX::Solver::Generic& s)
{
  bool ok;
  counter.incrementNumLineSearches();
  isNewtonDirection = false;
  NOX::Abstract::Vector& searchDirection = *tensorVecPtr;

  if ((counter.getNumLineSearches() == 1)  ||  (lsType == Newton))
  {
    isNewtonDirection = true;
    searchDirection = *newtonVecPtr;
  }

  // Do line search and compute new soln.
  if ((lsType != Dual) || (isNewtonDirection))
    ok = performLinesearch(newGrp, in_stepSize, searchDirection, s);
  else if (lsType == Dual)
  {
    double fTensor = 0.0;
    double fNew = 0.0;
    double tensorStep = 1.0;
    bool isTensorDescent = false;

    const Abstract::Group& oldGrp = s.getPreviousSolutionGroup();
    double fprime = slopeObj.computeSlope(searchDirection, oldGrp);

    // Backtrack along tensor direction if it is descent direction.
    if (fprime < 0)
    {
      ok = performLinesearch(newGrp, in_stepSize, searchDirection, s);
      assert(ok);
      fTensor = 0.5 * newGrp.getNormF() * newGrp.getNormF();
      tensorStep = in_stepSize;
      isTensorDescent = true;
    }

    // Backtrack along the Newton direction.
    ok = performLinesearch(newGrp, in_stepSize, *newtonVecPtr, s);
    fNew = 0.5 * newGrp.getNormF() * newGrp.getNormF();

    // If backtracking on the tensor step produced a better step, then use it.
    if (isTensorDescent  &&  (fTensor <= fNew))
    {
      newGrp.computeX(oldGrp, *tensorVecPtr, tensorStep);
      newGrp.computeF();
    }
  }

  return ok;
}
double NOX::MeritFunction::SumOfSquares::
computef(const NOX::Abstract::Group& grp) const
{
  if ( !(grp.isF()) ) {
    utils->err()
      << "ERROR: NOX::MeritFunction::SumOfSquares::computef() - "
      << "F has not been computed yet!.  Please call "
      << "computeF() on the group passed into this function."
      << std::endl;
    throw "NOX Error";
  }

  return (0.5 * grp.getNormF() * grp.getNormF());
}
Beispiel #3
0
double NOX::StatusTest::NormF::computeNorm(const NOX::Abstract::Group& grp)
{
  if (!grp.isF())
    return -1.0;

  double norm;
  int n = grp.getX().length();

  switch (normType) 
  {
    
  case NOX::Abstract::Vector::TwoNorm:
    norm = grp.getNormF();
    if (scaleType == Scaled)
      norm /= sqrt(1.0 * n);
    break;

  default:
    norm = grp.getF().norm(normType);
    if (scaleType == Scaled)
      norm /= n;
    break;

  }

  return norm;
}
bool  NOX::Direction::ModifiedNewton::rescueBadNewtonSolve(const NOX::Abstract::Group& grp) const
{
  //! Check if the "rescue" option has been selected
  if (!doRescue)
    return false;

  //! See if the group has compute the accuracy
  double accuracy;
  NOX::Abstract::Group::ReturnType status = oldJacobianGrpPtr->getNormLastLinearSolveResidual(accuracy);
    
  // If this functionality is not supported in the group, return false
  /* NOTE FROM TAMMY: We could later modify this to acutally caluclate
     the error itself if it's just a matter of the status being
     NotDefined. */
  if (status != NOX::Abstract::Group::Ok) 
    return false;

  // Check if there is any improvement in the relative residual
  double normF = grp.getNormF();

  // If we can't reduce the relative norm at all, we're not happy
  if (accuracy >= normF) 
    return false;

  // Otherwise, we just print a warning and keep going
  if (utils->isPrintType(NOX::Utils::Warning))
    utils->out() << "WARNING: NOX::Direction::ModifiedNewton::compute - Unable to achieve desired linear solve accuracy." << std::endl;
  return true;

}
void
NOX::Solver::TensorBased::printDirectionInfo(std::string dirName,
                    const NOX::Abstract::Vector& dir,
                    const NOX::Abstract::Group& soln,
                    bool isTensorModel) const
{
  double dirNorm = dir.norm();

  double residual = getNormModelResidual(dir, soln, isTensorModel);
  double residualRel = residual / soln.getNormF();

  double fprime = getDirectionalDerivative(dir, soln);
  double fprimeRel = fprime / dirNorm;

  if (utilsPtr->isPrintType(NOX::Utils::Details))
  {
    utilsPtr->out() << " " << dirName << " norm of model residual =   "
     << utilsPtr->sciformat(residual, 6) << " (abs)     "
     << utilsPtr->sciformat(residualRel, 6) << " (rel)" << std::endl;
    utilsPtr->out() << " " << dirName << " directional derivative =  "
     << utilsPtr->sciformat(fprime, 6) << " (abs)    "
     << utilsPtr->sciformat(fprimeRel, 6) << " (rel)" << std::endl;
    utilsPtr->out() << " " << dirName << " norm = "
       << utilsPtr->sciformat(dirNorm, 6) << std::endl;
  }
}
double NOX::LineSearch::Polynomial::
computeValue(const NOX::Abstract::Group& grp, double phi)
{
  double value = phi;

  if (suffDecrCond == AredPred)
  {
    value = grp.getNormF();
  }

  return value;
}
bool NOX::Direction::Broyden::doRestart(NOX::Abstract::Group& soln, 
					const NOX::Solver::LineSearchBased& solver)
{
  // Test 1 - First iteration!
  if (solver.getNumIterations() == 0)
    return true;

  // Test 2 - Frequency
  if (cnt >= cntMax)
    return true;

  // Test 3 - Last step was zero!
  if (solver.getStepSize() == 0.0)
    return true;

  // Test 4 - Check for convergence rate
  convRate = soln.getNormF() / solver.getPreviousSolutionGroup().getNormF();
  if (convRate > maxConvRate)
    return true;

  return false;
}
// **************************************************************************
// *** computeForcingTerm
// **************************************************************************
double NOX::Direction::Utils::InexactNewton::
computeForcingTerm(const NOX::Abstract::Group& soln,
		   const NOX::Abstract::Group& oldsoln, 
		   int niter,
		   const NOX::Solver::Generic& solver,
		   double eta_last)
{
  const std::string indent = "       ";

  if (forcingTermMethod == Constant) {
    if (printing->isPrintType(NOX::Utils::Details)) {
      printing->out() << indent << "CALCULATING FORCING TERM" << std::endl;
      printing->out() << indent << "Method: Constant" << std::endl;
      printing->out() << indent << "Forcing Term: " << eta_k << std::endl;
    }
    if (setTolerance)
      paramsPtr->sublist(directionMethod).sublist("Linear Solver").
	set("Tolerance", eta_k);

    return eta_k;
  }

  // Get linear solver current tolerance. 
  // NOTE: These values are changing at each nonlinear iteration and 
  // must either be updated from the parameter list each time a compute 
  // is called or supplied during the function call!
  double eta_km1 = 0.0;
  if (eta_last < 0.0)
    eta_km1 = paramsPtr->sublist(directionMethod).
      sublist("Linear Solver").get("Tolerance", 0.0);
  else
    eta_km1 = eta_last;

  // Tolerance may have been adjusted in a line search algorithm so we 
  // have to account for this.
  const NOX::Solver::LineSearchBased* solverPtr = 0;
  solverPtr = dynamic_cast<const NOX::Solver::LineSearchBased*>(&solver);
  if (solverPtr != 0) {
    eta_km1 = 1.0 - solverPtr->getStepSize() * (1.0 - eta_km1);
  }

  if (printing->isPrintType(NOX::Utils::Details)) {
    printing->out() << indent << "CALCULATING FORCING TERM" << std::endl;
    printing->out() << indent << "Method: " << method << std::endl;
  }


  if (forcingTermMethod == Type1) {
    
    if (niter == 0) {
      
      eta_k = eta_initial;

    }
    else {

      // Return norm of predicted F

      // do NOT use the following lines!! This does NOT account for 
      // line search step length taken.
      // const double normpredf = 0.0;
      // oldsoln.getNormLastLinearSolveResidual(normpredf);
      
      // Create a new vector to be the predicted RHS
      if (Teuchos::is_null(predRhs)) {
	predRhs = oldsoln.getF().clone(ShapeCopy);
      }
      if (Teuchos::is_null(stepDir)) {
	stepDir = oldsoln.getF().clone(ShapeCopy);
      }
      
      // stepDir = X - oldX (i.e., the step times the direction)
      stepDir->update(1.0, soln.getX(), -1.0, oldsoln.getX(), 0);
      
      // Compute predRhs = Jacobian * step * dir
      if (!(oldsoln.isJacobian())) {
	if (printing->isPrintType(NOX::Utils::Details)) {
	  printing->out() << "WARNING: NOX::InexactNewtonUtils::resetForcingTerm() - "
	       << "Jacobian is out of date! Recomputing Jacobian." << std::endl;
	}
	const_cast<NOX::Abstract::Group&>(oldsoln).computeJacobian();
      }
      oldsoln.applyJacobian(*stepDir, *predRhs);

      // Compute predRhs = RHSVector + predRhs (this is the predicted RHS)
      predRhs->update(1.0, oldsoln.getF(), 1.0);
      
      // Compute the norms
      double normpredf = predRhs->norm();
      double normf = soln.getNormF();
      double normoldf = oldsoln.getNormF();

      if (printing->isPrintType(NOX::Utils::Details)) {
	printing->out() << indent << "Forcing Term Norm: Using L-2 Norm."
			<< std::endl;
      }

      // Compute forcing term
      eta_k = fabs(normf - normpredf) / normoldf;
      
      // Some output
      if (printing->isPrintType(NOX::Utils::Details)) {
	printing->out() << indent << "Residual Norm k-1 =             " 
	     << normoldf << "\n";
	printing->out() << indent << "Residual Norm Linear Model k =  " 
	     << normpredf << "\n";
	printing->out() << indent << "Residual Norm k =               " << normf << "\n";
	printing->out() << indent << "Calculated eta_k (pre-bounds) = " << eta_k << std::endl;
      }
      
      // Impose safeguard and constraints ...
      const double tmp = (1.0 + sqrt(5.0)) / 2.0;
      const double eta_km1_alpha = pow(eta_km1, tmp);
      if (eta_km1_alpha > 0.1) 
	eta_k = NOX_MAX(eta_k, eta_km1_alpha);
      eta_k = NOX_MAX(eta_k, eta_min);
      eta_k = NOX_MIN(eta_max, eta_k);
    }
  }
    
  else if (forcingTermMethod == Type2) {  
    
    if (niter == 0) {
      
      eta_k = eta_initial;
      
    }
    else {

      double normf = soln.getNormF();
      double normoldf = oldsoln.getNormF();
      
      if (printing->isPrintType(NOX::Utils::Details)) {
	printing->out() << indent << "Forcing Term Norm: Using L-2 Norm."
			<< std::endl;
      }

      const double residual_ratio = normf / normoldf;
      
      eta_k = gamma * pow(residual_ratio, alpha);
      
      // Some output
      if (printing->isPrintType(NOX::Utils::Details)) {
	printing->out() << indent << "Residual Norm k-1 =             " 
			<< normoldf << "\n";
	printing->out() << indent << "Residual Norm k =               " 
			<< normf << "\n";
	printing->out() << indent << "Calculated eta_k (pre-bounds) = " 
			<< eta_k << std::endl;
      }
      
      // Impose safeguard and constraints ... 
      const double eta_k_alpha = gamma * pow(eta_km1, alpha);
      if (eta_k_alpha > 0.1) 
	eta_k = NOX_MAX(eta_k, eta_k_alpha);
      eta_k = NOX_MAX(eta_k, eta_min);
      eta_k = NOX_MIN(eta_max, eta_k);
    }
    
  }
  
  // Set the new linear solver tolerance
  if (setTolerance) 
    paramsPtr->sublist(directionMethod).sublist("Linear Solver").
      set("Tolerance", eta_k);

  if (printing->isPrintType(NOX::Utils::Details)) 
    printing->out() << indent << "Forcing Term: " << eta_k << std::endl;
  
  return eta_k;
}
bool
NOX::Solver::TensorBased::performLinesearch(NOX::Abstract::Group& newSoln,
                        double& in_stepSize,
                        const NOX::Abstract::Vector& lsDir,
                        const NOX::Solver::Generic& s)
{
  if (print.isPrintType(NOX::Utils::InnerIteration))
  {
    utilsPtr->out() << "\n" << NOX::Utils::fill(72) << "\n";
    utilsPtr->out() << "-- Tensor Line Search (";
    if (lsType == Curvilinear)
      utilsPtr->out() << "Curvilinear";
    else if (lsType == Standard)
      utilsPtr->out() << "Standard";
    else if (lsType == FullStep)
      utilsPtr->out() << "Full Step";
    else if (lsType == Dual)
      utilsPtr->out() << "Dual";
    utilsPtr->out() << ") -- " << std::endl;
  }

  // Local variables
  bool isFailed = false;
  bool isAcceptable = false;
  bool isFirstPass = true;
  std::string message = "(STEP ACCEPTED!)";

  // Set counters
  int lsIterations = 1;

  // Get Old f
  const Abstract::Group& oldSoln = s.getPreviousSolutionGroup();
  double fOld = 0.5 * oldSoln.getNormF() * oldSoln.getNormF();

  // Compute first trial point and its function value
  in_stepSize = defaultStep;
  newSoln.computeX(oldSoln, lsDir, in_stepSize);
  newSoln.computeF();
  double fNew = 0.5 * newSoln.getNormF() * newSoln.getNormF();

  // Stop here if only using the full step
  if (lsType == FullStep)
  {
      print.printStep(lsIterations, in_stepSize, fOld, fNew, message);
      return (!isFailed);
  }

  // Compute directional derivative
  double fprime;
  if ((lsType == Curvilinear)  &&  !(isNewtonDirection))
    fprime = slopeObj.computeSlope(*newtonVecPtr, oldSoln);
  else
    fprime = slopeObj.computeSlope(lsDir, oldSoln);
  numJvMults++;  // computeSlope() has J*v inside of it

  // Compute the convergence criteria for the line search
  double threshold = fOld + alpha*in_stepSize*fprime;
  isAcceptable = (fNew < threshold);

  // Update counter and temporarily hold direction if a linesearch is needed
  if (!isAcceptable)
  {
    counter.incrementNumNonTrivialLineSearches();
    *tmpVecPtr = lsDir;
  }

  // Iterate until the trial point is accepted....
  while (!isAcceptable)
  {
    // Check for linesearch failure
    if (lsIterations > maxIters)
    {
      isFailed = true;
      message = "(FAILED - Max Iters)";
      break;
    }

    print.printStep(lsIterations, in_stepSize, fOld, fNew);

    // Is the full tensor step a descent direction?  If not, switch to Newton
    if (isFirstPass &&
    (!isNewtonDirection) &&
    (fprime >= 0) &&
    (lsType != Curvilinear) )
    {
      *tmpVecPtr = *newtonVecPtr;
      fprime = slopeObj.computeSlope(*tmpVecPtr, oldSoln);
      numJvMults++;

      if (utilsPtr->isPrintType(NOX::Utils::Details))
    utilsPtr->out() << "  Switching to Newton step.  New fprime = "
         << utilsPtr->sciformat(fprime, 6) << std::endl;
    }
    else
    {
      in_stepSize = selectLambda(fNew, fOld, fprime, in_stepSize);
    }

    isFirstPass = false;

    // Check for linesearch failure
    if (in_stepSize < minStep)
    {
      isFailed = true;
      message = "(FAILED - Min Step)";
      break;
    }

    // Update the number of linesearch iterations
    counter.incrementNumIterations();
    lsIterations ++;

    // Compute new trial point and its function value
    if ((lsType == Curvilinear) && !(isNewtonDirection))
    {
      computeCurvilinearStep(*tmpVecPtr, oldSoln, s, in_stepSize);
      // Note: oldSoln is needed above to get correct preconditioner
      newSoln.computeX(oldSoln, *tmpVecPtr, 1.0);
    }
    else
    {
      newSoln.computeX(oldSoln, *tmpVecPtr, in_stepSize);
    }
    newSoln.computeF();
    fNew = 0.5 * newSoln.getNormF() * newSoln.getNormF();

    // Recompute convergence criteria based on new step
    threshold = fOld + alpha*in_stepSize*fprime;
    isAcceptable = (fNew < threshold);
  }


  if (isFailed)
  {
    counter.incrementNumFailedLineSearches();

    if (recoveryStepType == Constant)
    {
      in_stepSize = recoveryStep;
      if (in_stepSize == 0.0)
      {
    newSoln = oldSoln;
    newSoln.computeF();
    fNew = fOld;
      }
      else
      {
    // Update the group using recovery step
    if ((lsType == Curvilinear) && !(isNewtonDirection))
    {
      computeCurvilinearStep(*tmpVecPtr, oldSoln, s, in_stepSize);
      // Note: oldSoln is needed above to get correct preconditioner
      newSoln.computeX(oldSoln, *tmpVecPtr, 1.0);
    }
    else
    {
      newSoln.computeX(oldSoln, *tmpVecPtr, in_stepSize);
    }
    //newSoln.computeX(oldSoln, lsDir, in_stepSize);
    newSoln.computeF();
    fNew = 0.5 * newSoln.getNormF() * newSoln.getNormF();
    message = "(USING RECOVERY STEP!)";
      }
    }
    else
      message = "(USING LAST STEP!)";
  }

  print.printStep(lsIterations, in_stepSize, fOld, fNew, message);
  counter.setValues(paramsPtr->sublist("Line Search"));

  return (!isFailed);
}