Number PDFullSpaceSolver::ComputeResidualRatio(const IteratesVector& rhs,
      const IteratesVector& res,
      const IteratesVector& resid)
  {
    DBG_START_METH("PDFullSpaceSolver::ComputeResidualRatio", dbg_verbosity);

    Number nrm_rhs = rhs.Amax();
    Number nrm_res = res.Amax();
    Number nrm_resid = resid.Amax();
    Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                   "nrm_rhs = %8.2e nrm_sol = %8.2e nrm_resid = %8.2e\n",
                   nrm_rhs, nrm_res, nrm_resid);

    if (nrm_rhs+nrm_res == 0.) {
      return nrm_resid;  // this should be zero
    }
    else {
      // ToDo: determine how to include norm of matrix, and what
      // safeguard to use against incredibly large solution vectors
      Number max_cond = 1e6;
      return nrm_resid/(Min(nrm_res, max_cond*nrm_rhs)+nrm_rhs);
    }
  }
  void PDFullSpaceSolver::ComputeResiduals(
    const SymMatrix& W,
    const Matrix& J_c,
    const Matrix& J_d,
    const Matrix& Px_L,
    const Matrix& Px_U,
    const Matrix& Pd_L,
    const Matrix& Pd_U,
    const Vector& z_L,
    const Vector& z_U,
    const Vector& v_L,
    const Vector& v_U,
    const Vector& slack_x_L,
    const Vector& slack_x_U,
    const Vector& slack_s_L,
    const Vector& slack_s_U,
    const Vector& sigma_x,
    const Vector& sigma_s,
    Number alpha,
    Number beta,
    const IteratesVector& rhs,
    const IteratesVector& res,
    IteratesVector& resid)
  {
    DBG_START_METH("PDFullSpaceSolver::ComputeResiduals", dbg_verbosity);

    DBG_PRINT_VECTOR(2, "res", res);
    IpData().TimingStats().ComputeResiduals().Start();

    // Get the current sizes of the perturbation factors
    Number delta_x;
    Number delta_s;
    Number delta_c;
    Number delta_d;
    perturbHandler_->CurrentPerturbation(delta_x, delta_s, delta_c, delta_d);

    SmartPtr<Vector> tmp;

    // x
    W.MultVector(1., *res.x(), 0., *resid.x_NonConst());
    J_c.TransMultVector(1., *res.y_c(), 1., *resid.x_NonConst());
    J_d.TransMultVector(1., *res.y_d(), 1., *resid.x_NonConst());
    Px_L.MultVector(-1., *res.z_L(), 1., *resid.x_NonConst());
    Px_U.MultVector(1., *res.z_U(), 1., *resid.x_NonConst());
    resid.x_NonConst()->AddTwoVectors(delta_x, *res.x(), -1., *rhs.x(), 1.);

    // s
    Pd_U.MultVector(1., *res.v_U(), 0., *resid.s_NonConst());
    Pd_L.MultVector(-1., *res.v_L(), 1., *resid.s_NonConst());
    resid.s_NonConst()->AddTwoVectors(-1., *res.y_d(), -1., *rhs.s(), 1.);
    if (delta_s!=0.) {
      resid.s_NonConst()->Axpy(delta_s, *res.s());
    }

    // c
    J_c.MultVector(1., *res.x(), 0., *resid.y_c_NonConst());
    resid.y_c_NonConst()->AddTwoVectors(-delta_c, *res.y_c(), -1., *rhs.y_c(), 1.);

    // d
    J_d.MultVector(1., *res.x(), 0., *resid.y_d_NonConst());
    resid.y_d_NonConst()->AddTwoVectors(-1., *res.s(), -1., *rhs.y_d(), 1.);
    if (delta_d!=0.) {
      resid.y_d_NonConst()->Axpy(-delta_d, *res.y_d());
    }

    // zL
    resid.z_L_NonConst()->Copy(*res.z_L());
    resid.z_L_NonConst()->ElementWiseMultiply(slack_x_L);
    tmp = z_L.MakeNew();
    Px_L.TransMultVector(1., *res.x(), 0., *tmp);
    tmp->ElementWiseMultiply(z_L);
    resid.z_L_NonConst()->AddTwoVectors(1., *tmp, -1., *rhs.z_L(), 1.);

    // zU
    resid.z_U_NonConst()->Copy(*res.z_U());
    resid.z_U_NonConst()->ElementWiseMultiply(slack_x_U);
    tmp = z_U.MakeNew();
    Px_U.TransMultVector(1., *res.x(), 0., *tmp);
    tmp->ElementWiseMultiply(z_U);
    resid.z_U_NonConst()->AddTwoVectors(-1., *tmp, -1., *rhs.z_U(), 1.);

    // vL
    resid.v_L_NonConst()->Copy(*res.v_L());
    resid.v_L_NonConst()->ElementWiseMultiply(slack_s_L);
    tmp = v_L.MakeNew();
    Pd_L.TransMultVector(1., *res.s(), 0., *tmp);
    tmp->ElementWiseMultiply(v_L);
    resid.v_L_NonConst()->AddTwoVectors(1., *tmp, -1., *rhs.v_L(), 1.);

    // vU
    resid.v_U_NonConst()->Copy(*res.v_U());
    resid.v_U_NonConst()->ElementWiseMultiply(slack_s_U);
    tmp = v_U.MakeNew();
    Pd_U.TransMultVector(1., *res.s(), 0., *tmp);
    tmp->ElementWiseMultiply(v_U);
    resid.v_U_NonConst()->AddTwoVectors(-1., *tmp, -1., *rhs.v_U(), 1.);

    DBG_PRINT_VECTOR(2, "resid", resid);

    if (Jnlst().ProduceOutput(J_MOREVECTOR, J_LINEAR_ALGEBRA)) {
      resid.Print(Jnlst(), J_MOREVECTOR, J_LINEAR_ALGEBRA, "resid");
    }

    if (Jnlst().ProduceOutput(J_MOREDETAILED, J_LINEAR_ALGEBRA)) {
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_x  %e\n", resid.x()->Amax());
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_s  %e\n", resid.s()->Amax());
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_c  %e\n", resid.y_c()->Amax());
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_d  %e\n", resid.y_d()->Amax());
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_zL %e\n", resid.z_L()->Amax());
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_zU %e\n", resid.z_U()->Amax());
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_vL %e\n", resid.v_L()->Amax());
      Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                     "max-norm resid_vU %e\n", resid.v_U()->Amax());
    }
    IpData().TimingStats().ComputeResiduals().End();
  }
  bool PDFullSpaceSolver::SolveOnce(bool resolve_with_better_quality,
                                    bool pretend_singular,
                                    const SymMatrix& W,
                                    const Matrix& J_c,
                                    const Matrix& J_d,
                                    const Matrix& Px_L,
                                    const Matrix& Px_U,
                                    const Matrix& Pd_L,
                                    const Matrix& Pd_U,
                                    const Vector& z_L,
                                    const Vector& z_U,
                                    const Vector& v_L,
                                    const Vector& v_U,
                                    const Vector& slack_x_L,
                                    const Vector& slack_x_U,
                                    const Vector& slack_s_L,
                                    const Vector& slack_s_U,
                                    const Vector& sigma_x,
                                    const Vector& sigma_s,
                                    Number alpha,
                                    Number beta,
                                    const IteratesVector& rhs,
                                    IteratesVector& res)
  {
    // TO DO LIST:
    //
    // 1. decide for reasonable return codes (e.g. fatal error, too
    //    ill-conditioned...)
    // 2. Make constants parameters that can be set from the outside
    // 3. Get Information out of Ipopt structures
    // 4. add heuristic for structurally singular problems
    // 5. see if it makes sense to distinguish delta_x and delta_s,
    //    or delta_c and delta_d
    // 6. increase pivot tolerance if number of get evals so too small
    DBG_START_METH("PDFullSpaceSolver::SolveOnce",dbg_verbosity);

    IpData().TimingStats().PDSystemSolverSolveOnce().Start();

    // Compute the right hand side for the augmented system formulation
    SmartPtr<Vector> augRhs_x = rhs.x()->MakeNewCopy();
    Px_L.AddMSinvZ(1.0, slack_x_L, *rhs.z_L(), *augRhs_x);
    Px_U.AddMSinvZ(-1.0, slack_x_U, *rhs.z_U(), *augRhs_x);

    SmartPtr<Vector> augRhs_s = rhs.s()->MakeNewCopy();
    Pd_L.AddMSinvZ(1.0, slack_s_L, *rhs.v_L(), *augRhs_s);
    Pd_U.AddMSinvZ(-1.0, slack_s_U, *rhs.v_U(), *augRhs_s);

    // Get space into which we can put the solution of the augmented system
    SmartPtr<IteratesVector> sol = res.MakeNewIteratesVector(true);

    // Now check whether any data has changed
    std::vector<const TaggedObject*> deps(13);
    deps[0] = &W;
    deps[1] = &J_c;
    deps[2] = &J_d;
    deps[3] = &z_L;
    deps[4] = &z_U;
    deps[5] = &v_L;
    deps[6] = &v_U;
    deps[7] = &slack_x_L;
    deps[8] = &slack_x_U;
    deps[9] = &slack_s_L;
    deps[10] = &slack_s_U;
    deps[11] = &sigma_x;
    deps[12] = &sigma_s;
    void* dummy;
    bool uptodate = dummy_cache_.GetCachedResult(dummy, deps);
    if (!uptodate) {
      dummy_cache_.AddCachedResult(dummy, deps);
      augsys_improved_ = false;
    }
    // improve_current_solution can only be true, if that system has
    // been solved before
    DBG_ASSERT((!resolve_with_better_quality && !pretend_singular) || uptodate);

    ESymSolverStatus retval;
    if (uptodate && !pretend_singular) {

      // Get the perturbation values
      Number delta_x;
      Number delta_s;
      Number delta_c;
      Number delta_d;
      perturbHandler_->CurrentPerturbation(delta_x, delta_s, delta_c, delta_d);

      // No need to go through the pain of finding the appropriate
      // values for the deltas, because the matrix hasn't changed since
      // the last call.  So, just call the Solve Method
      //
      // Note: resolve_with_better_quality is true, then the Solve
      // method has already asked the augSysSolver to increase the
      // quality at the end solve, and we are now getting the solution
      // with that better quality
      retval = augSysSolver_->Solve(&W, 1.0, &sigma_x, delta_x,
                                    &sigma_s, delta_s, &J_c, NULL,
                                    delta_c, &J_d, NULL, delta_d,
                                    *augRhs_x, *augRhs_s, *rhs.y_c(), *rhs.y_d(),
                                    *sol->x_NonConst(), *sol->s_NonConst(),
                                    *sol->y_c_NonConst(), *sol->y_d_NonConst(),
                                    false, 0);
      if (retval!=SYMSOLVER_SUCCESS) {
        IpData().TimingStats().PDSystemSolverSolveOnce().End();
        return false;
      }
    }
    else {
      const Index numberOfEVals=rhs.y_c()->Dim()+rhs.y_d()->Dim();
      // counter for the number of trial evaluations
      // (ToDo is not at the correct place)
      Index count = 0;

      // Get the very first perturbation values from the perturbation
      // Handler
      Number delta_x;
      Number delta_s;
      Number delta_c;
      Number delta_d;
      perturbHandler_->ConsiderNewSystem(delta_x, delta_s, delta_c, delta_d);

      retval = SYMSOLVER_SINGULAR;
      bool fail = false;

      while (retval!= SYMSOLVER_SUCCESS && !fail) {

        if (pretend_singular) {
          retval = SYMSOLVER_SINGULAR;
          pretend_singular = false;
        }
        else {
          count++;
          Jnlst().Printf(J_MOREDETAILED, J_LINEAR_ALGEBRA,
                         "Solving system with delta_x=%e delta_s=%e\n                    delta_c=%e delta_d=%e\n",
                         delta_x, delta_s, delta_c, delta_d);
          bool check_inertia = true;
          if (neg_curv_test_tol_ > 0.) {
            check_inertia = false;
          }
          retval = augSysSolver_->Solve(&W, 1.0, &sigma_x, delta_x,
                                        &sigma_s, delta_s, &J_c, NULL,
                                        delta_c, &J_d, NULL, delta_d,
                                        *augRhs_x, *augRhs_s, *rhs.y_c(), *rhs.y_d(),
                                        *sol->x_NonConst(), *sol->s_NonConst(),
                                        *sol->y_c_NonConst(), *sol->y_d_NonConst(),                                     check_inertia, numberOfEVals);
        }
        if (retval==SYMSOLVER_FATAL_ERROR) return false;
        if (retval==SYMSOLVER_SINGULAR &&
            (rhs.y_c()->Dim()+rhs.y_d()->Dim() > 0) ) {

          // Get new perturbation factors from the perturbation
          // handlers for the singular case
          bool pert_return = perturbHandler_->PerturbForSingularity(delta_x, delta_s,
                             delta_c, delta_d);
          if (!pert_return) {
            Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                           "PerturbForSingularity can't be done\n");
            IpData().TimingStats().PDSystemSolverSolveOnce().End();
            return false;
          }
        }
        else if (retval==SYMSOLVER_WRONG_INERTIA &&
                 augSysSolver_->NumberOfNegEVals() < numberOfEVals) {
          Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                         "Number of negative eigenvalues too small!\n");
          // If the number of negative eigenvalues is too small, then
          // we first try to remedy this by asking for better quality
          // solution (e.g. increasing pivot tolerance), and if that
          // doesn't help, we assume that the system is singular
          bool assume_singular = true;
          if (!augsys_improved_) {
            Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                           "Asking augmented system solver to improve quality of its solutions.\n");
            augsys_improved_ = augSysSolver_->IncreaseQuality();
            if (augsys_improved_) {
              IpData().Append_info_string("q");
              assume_singular = false;
            }
            else {
              Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                             "Quality could not be improved\n");
            }
          }
          if (assume_singular) {
            bool pert_return =
                               perturbHandler_->PerturbForSingularity(delta_x, delta_s,
                                                                      delta_c, delta_d);
            if (!pert_return) {
              Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                             "PerturbForSingularity can't be done for assume singular.\n");
              IpData().TimingStats().PDSystemSolverSolveOnce().End();
              return false;
            }
            IpData().Append_info_string("a");
          }
        }
        else if (retval==SYMSOLVER_WRONG_INERTIA ||
                 retval==SYMSOLVER_SINGULAR) {
          // Get new perturbation factors from the perturbation
          // handlers for the case of wrong inertia
          bool pert_return = perturbHandler_->PerturbForWrongInertia(delta_x, delta_s,
                             delta_c, delta_d);
          if (!pert_return) {
            Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                           "PerturbForWrongInertia can't be done for wrong interia or singular.\n");
            IpData().TimingStats().PDSystemSolverSolveOnce().End();
            return false;
          }
        }
        else if (neg_curv_test_tol_ > 0.) {
          DBG_ASSERT(augSysSolver_->ProvidesInertia());
          // we now check if the inertia is possible wrong
          Index neg_values = augSysSolver_->NumberOfNegEVals();
          if (neg_values != numberOfEVals) {
            // check if we have a direction of sufficient positive curvature
            SmartPtr<Vector> x_tmp = sol->x()->MakeNew();
            W.MultVector(1., *sol->x(), 0., *x_tmp);
            Number xWx = x_tmp->Dot(*sol->x());
            x_tmp->Copy(*sol->x());
            x_tmp->ElementWiseMultiply(sigma_x);
            xWx += x_tmp->Dot(*sol->x());
            SmartPtr<Vector> s_tmp = sol->s()->MakeNewCopy();
            s_tmp->ElementWiseMultiply(sigma_s);
            xWx += s_tmp->Dot(*sol->s());
            if (neg_curv_test_reg_) {
              x_tmp->Copy(*sol->x());
              x_tmp->Scal(delta_x);
              xWx += x_tmp->Dot(*sol->x());

              s_tmp->Copy(*sol->s());
              s_tmp->Scal(delta_s);
              xWx += s_tmp->Dot(*sol->s());
            }
            Number xs_nrmsq = pow(sol->x()->Nrm2(),2) + pow(sol->s()->Nrm2(),2);
            Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                           "In inertia heuristic: xWx = %e xx = %e\n",
                           xWx, xs_nrmsq);
            if (xWx < neg_curv_test_tol_*xs_nrmsq) {
              Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                             "    -> Redo with modified matrix.\n");
              bool pert_return = perturbHandler_->PerturbForWrongInertia(delta_x, delta_s,
                                 delta_c, delta_d);
              if (!pert_return) {
                Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                               "PerturbForWrongInertia can't be done for inertia heuristic.\n");
                IpData().TimingStats().PDSystemSolverSolveOnce().End();
                return false;
              }
              retval = SYMSOLVER_WRONG_INERTIA;
            }
          }
        }
      } // while (retval!=SYMSOLVER_SUCCESS && !fail) {

      // Some output
      Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                     "Number of trial factorizations performed: %d\n",
                     count);
      Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                     "Perturbation parameters: delta_x=%e delta_s=%e\n                         delta_c=%e delta_d=%e\n",
                     delta_x, delta_s, delta_c, delta_d);
      // Set the perturbation values in the Data object
      IpData().setPDPert(delta_x, delta_s, delta_c, delta_d);
    }

    // Compute the remaining sol Vectors
    Px_L.SinvBlrmZMTdBr(-1., slack_x_L, *rhs.z_L(), z_L, *sol->x(), *sol->z_L_NonConst());
    Px_U.SinvBlrmZMTdBr(1., slack_x_U, *rhs.z_U(), z_U, *sol->x(), *sol->z_U_NonConst());
    Pd_L.SinvBlrmZMTdBr(-1., slack_s_L, *rhs.v_L(), v_L, *sol->s(), *sol->v_L_NonConst());
    Pd_U.SinvBlrmZMTdBr(1., slack_s_U, *rhs.v_U(), v_U, *sol->s(), *sol->v_U_NonConst());

    // Finally let's assemble the res result vectors
    res.AddOneVector(alpha, *sol, beta);

    IpData().TimingStats().PDSystemSolverSolveOnce().End();

    return true;
  }
  bool PDFullSpaceSolver::Solve(Number alpha,
                                Number beta,
                                const IteratesVector& rhs,
                                IteratesVector& res,
                                bool allow_inexact,
                                bool improve_solution /* = false */)
  {
    DBG_START_METH("PDFullSpaceSolver::Solve",dbg_verbosity);
    DBG_ASSERT(!allow_inexact || !improve_solution);
    DBG_ASSERT(!improve_solution || beta==0.);

    // Timing of PDSystem solver starts here
    IpData().TimingStats().PDSystemSolverTotal().Start();

    DBG_PRINT_VECTOR(2, "rhs_x", *rhs.x());
    DBG_PRINT_VECTOR(2, "rhs_s", *rhs.s());
    DBG_PRINT_VECTOR(2, "rhs_c", *rhs.y_c());
    DBG_PRINT_VECTOR(2, "rhs_d", *rhs.y_d());
    DBG_PRINT_VECTOR(2, "rhs_zL", *rhs.z_L());
    DBG_PRINT_VECTOR(2, "rhs_zU", *rhs.z_U());
    DBG_PRINT_VECTOR(2, "rhs_vL", *rhs.v_L());
    DBG_PRINT_VECTOR(2, "rhs_vU", *rhs.v_U());
    DBG_PRINT_VECTOR(2, "res_x in", *res.x());
    DBG_PRINT_VECTOR(2, "res_s in", *res.s());
    DBG_PRINT_VECTOR(2, "res_c in", *res.y_c());
    DBG_PRINT_VECTOR(2, "res_d in", *res.y_d());
    DBG_PRINT_VECTOR(2, "res_zL in", *res.z_L());
    DBG_PRINT_VECTOR(2, "res_zU in", *res.z_U());
    DBG_PRINT_VECTOR(2, "res_vL in", *res.v_L());
    DBG_PRINT_VECTOR(2, "res_vU in", *res.v_U());

    // if beta is nonzero, keep a copy of the incoming values in res_ */
    SmartPtr<IteratesVector> copy_res;
    if (beta != 0.) {
      copy_res = res.MakeNewIteratesVectorCopy();
    }

    // Receive data about matrix
    SmartPtr<const Vector> x = IpData().curr()->x();
    SmartPtr<const Vector> s = IpData().curr()->s();
    SmartPtr<const SymMatrix> W = IpData().W();
    SmartPtr<const Matrix> J_c = IpCq().curr_jac_c();
    SmartPtr<const Matrix> J_d = IpCq().curr_jac_d();
    SmartPtr<const Matrix> Px_L = IpNLP().Px_L();
    SmartPtr<const Matrix> Px_U = IpNLP().Px_U();
    SmartPtr<const Matrix> Pd_L = IpNLP().Pd_L();
    SmartPtr<const Matrix> Pd_U = IpNLP().Pd_U();
    SmartPtr<const Vector> z_L = IpData().curr()->z_L();
    SmartPtr<const Vector> z_U = IpData().curr()->z_U();
    SmartPtr<const Vector> v_L = IpData().curr()->v_L();
    SmartPtr<const Vector> v_U = IpData().curr()->v_U();
    SmartPtr<const Vector> slack_x_L = IpCq().curr_slack_x_L();
    SmartPtr<const Vector> slack_x_U = IpCq().curr_slack_x_U();
    SmartPtr<const Vector> slack_s_L = IpCq().curr_slack_s_L();
    SmartPtr<const Vector> slack_s_U = IpCq().curr_slack_s_U();
    SmartPtr<const Vector> sigma_x = IpCq().curr_sigma_x();
    SmartPtr<const Vector> sigma_s = IpCq().curr_sigma_s();
    DBG_PRINT_VECTOR(2, "Sigma_x", *sigma_x);
    DBG_PRINT_VECTOR(2, "Sigma_s", *sigma_s);

    bool done = false;
    // The following flag is set to true, if we asked the linear
    // solver to improve the quality of the solution in
    // the next solve
    bool resolve_with_better_quality = false;
    // the following flag is set to true, if iterative refinement
    // failed and we want to try if a modified system is able to
    // remedy that problem by pretending the matrix is singular
    bool pretend_singular = false;
    bool pretend_singular_last_time = false;

    // Beginning of loop for solving the system (including all
    // modifications for the linear system to ensure good solution
    // quality)
    while (!done) {

      // if improve_solution is true, we are given already a solution
      // from the calling function, so we can skip the first solve
      bool solve_retval = true;
      if (!improve_solution) {
        solve_retval =
          SolveOnce(resolve_with_better_quality, pretend_singular,
                    *W, *J_c, *J_d, *Px_L, *Px_U, *Pd_L, *Pd_U, *z_L, *z_U,
                    *v_L, *v_U, *slack_x_L, *slack_x_U, *slack_s_L, *slack_s_U,
                    *sigma_x, *sigma_s, 1., 0., rhs, res);
        resolve_with_better_quality = false;
        pretend_singular = false;
      }
      improve_solution = false;

      if (!solve_retval) {
        // If system seems not to be solvable, we return with false
        // and let the calling routine deal with it.
        IpData().TimingStats().PDSystemSolverTotal().End();
        return false;
      }

      if (allow_inexact) {
        // no safety checks required
        if (Jnlst().ProduceOutput(J_MOREDETAILED, J_LINEAR_ALGEBRA)) {
          SmartPtr<IteratesVector> resid = res.MakeNewIteratesVector(true);
          ComputeResiduals(*W, *J_c, *J_d, *Px_L, *Px_U, *Pd_L, *Pd_U,
                           *z_L, *z_U, *v_L, *v_U, *slack_x_L, *slack_x_U,
                           *slack_s_L, *slack_s_U, *sigma_x, *sigma_s,
                           alpha, beta, rhs, res, *resid);
        }
        break;
      }

      // Get space for the residual
      SmartPtr<IteratesVector> resid = res.MakeNewIteratesVector(true);

      // ToDo don't to that after max refinement?
      ComputeResiduals(*W, *J_c, *J_d, *Px_L, *Px_U, *Pd_L, *Pd_U,
                       *z_L, *z_U, *v_L, *v_U, *slack_x_L, *slack_x_U,
                       *slack_s_L, *slack_s_U, *sigma_x, *sigma_s,
                       alpha, beta, rhs, res, *resid);

      Number residual_ratio =
        ComputeResidualRatio(rhs, res, *resid);
      Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                     "residual_ratio = %e\n", residual_ratio);
      Number residual_ratio_old = residual_ratio;

      // Beginning of loop for iterative refinement
      Index num_iter_ref = 0;
      bool quit_refinement = false;
      while (!allow_inexact && !quit_refinement &&
             (num_iter_ref < min_refinement_steps_ ||
              residual_ratio > residual_ratio_max_) ) {

        // To the next back solve
        solve_retval =
          SolveOnce(resolve_with_better_quality, false,
                    *W, *J_c, *J_d, *Px_L, *Px_U, *Pd_L, *Pd_U, *z_L, *z_U,
                    *v_L, *v_U, *slack_x_L, *slack_x_U, *slack_s_L, *slack_s_U,
                    *sigma_x, *sigma_s, -1., 1., *resid, res);
        ASSERT_EXCEPTION(solve_retval, INTERNAL_ABORT,
                         "SolveOnce returns false during iterative refinement.");

        ComputeResiduals(*W, *J_c, *J_d, *Px_L, *Px_U, *Pd_L, *Pd_U,
                         *z_L, *z_U, *v_L, *v_U, *slack_x_L, *slack_x_U,
                         *slack_s_L, *slack_s_U, *sigma_x, *sigma_s,
                         alpha, beta, rhs, res, *resid);

        residual_ratio =
          ComputeResidualRatio(rhs, res, *resid);
        Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                       "residual_ratio = %e\n", residual_ratio);

        num_iter_ref++;
        // Check if we have to give up on iterative refinement
        if (residual_ratio > residual_ratio_max_ &&
            num_iter_ref>min_refinement_steps_ &&
            (num_iter_ref>max_refinement_steps_ ||
             residual_ratio>residual_improvement_factor_*residual_ratio_old)) {

          Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                         "Iterative refinement failed with residual_ratio = %e\n", residual_ratio);
          quit_refinement = true;

          // Pretend singularity only once - if it didn't help, we
          // have to live with what we got so far
          resolve_with_better_quality = false;
          DBG_PRINT((1, "pretend_singular = %d\n", pretend_singular));
          if (!pretend_singular_last_time) {
            // First try if we can ask the augmented system solver to
            // improve the quality of the solution (only if that hasn't
            // been done before for this linear system)
            if (!augsys_improved_) {
              Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                             "Asking augmented system solver to improve quality of its solutions.\n");
              augsys_improved_ = augSysSolver_->IncreaseQuality();
              if (augsys_improved_) {
                IpData().Append_info_string("q");
                resolve_with_better_quality = true;
              }
              else {
                // solver said it cannot improve quality, so let
                // possibly conclude that the current modification is
                // singular
                pretend_singular = true;
              }
            }
            else {
              // we had already asked the solver before to improve the
              // quality of the solution, so let's now pretend that the
              // modification is possibly singular
              pretend_singular = true;
            }
            pretend_singular_last_time = pretend_singular;
            if (pretend_singular) {
              // let's only conclude that the current linear system
              // including modifications is singular, if the residual is
              // quite bad
              if (residual_ratio < residual_ratio_singular_) {
                pretend_singular = false;
                IpData().Append_info_string("S");
                Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                               "Just accept current solution.\n");
              }
              else {
                IpData().Append_info_string("s");
                Jnlst().Printf(J_DETAILED, J_LINEAR_ALGEBRA,
                               "Pretend that the current system (including modifications) is singular.\n");
              }
            }
          }
          else {
            pretend_singular = false;
            DBG_PRINT((1,"Resetting pretend_singular to false.\n"));
          }
        }

        residual_ratio_old = residual_ratio;
      } // End of loop for iterative refinement

      done = !(resolve_with_better_quality) && !(pretend_singular);

    } // End of loop for solving the linear system (incl. modifications)

    // Finally let's assemble the res result vectors
    if (alpha != 0.) {
      res.Scal(alpha);
    }

    if (beta != 0.) {
      res.Axpy(beta, *copy_res);
    }

    DBG_PRINT_VECTOR(2, "res_x", *res.x());
    DBG_PRINT_VECTOR(2, "res_s", *res.s());
    DBG_PRINT_VECTOR(2, "res_c", *res.y_c());
    DBG_PRINT_VECTOR(2, "res_d", *res.y_d());
    DBG_PRINT_VECTOR(2, "res_zL", *res.z_L());
    DBG_PRINT_VECTOR(2, "res_zU", *res.z_U());
    DBG_PRINT_VECTOR(2, "res_vL", *res.v_L());
    DBG_PRINT_VECTOR(2, "res_vU", *res.v_U());

    IpData().TimingStats().PDSystemSolverTotal().End();

    return true;
  }
  bool StdStepCalculator::Step(DenseVector& delta_u,
			       IteratesVector& sol)
  {
    DBG_START_METH("StdStepCalculator::Step", dbg_verbosity);

    bool retval = false;

    SmartPtr<IteratesVector> delta_u_long = IpData().trial()->MakeNewIteratesVector();
    ift_data_->TransMultiply(delta_u, *delta_u_long);

    SmartPtr<IteratesVector> r_s = IpData().trial()->MakeNewIteratesVector();
    if (kkt_residuals_) {
      /* This should be almost zero... */
      r_s->Set_x_NonConst(*IpCq().curr_grad_lag_x()->MakeNewCopy());
      r_s->Set_s_NonConst(*IpCq().curr_grad_lag_s()->MakeNewCopy());
      r_s->Set_y_c_NonConst(*IpCq().curr_c()->MakeNewCopy());
      r_s->Set_y_d_NonConst(*IpCq().curr_d_minus_s()->MakeNewCopy());
      r_s->Set_z_L_NonConst(*IpCq().curr_compl_x_L()->MakeNewCopy());
      r_s->Set_z_U_NonConst(*IpCq().curr_compl_x_U()->MakeNewCopy());
      r_s->Set_v_L_NonConst(*IpCq().curr_compl_s_L()->MakeNewCopy());
      r_s->Set_v_U_NonConst(*IpCq().curr_compl_s_U()->MakeNewCopy());

      r_s->Print(Jnlst(),J_VECTOR,J_USER1,"r_s init");
      delta_u.Print(Jnlst(),J_VECTOR,J_USER1,"delta_u init");
      DBG_PRINT((dbg_verbosity,"r_s init Nrm2=%23.16e\n", r_s->Asum()));

      delta_u_long->Axpy(-1.0, *r_s);
    }

    backsolver_->Solve(&sol, ConstPtr(delta_u_long));

    SmartPtr<IteratesVector> Kr_s;
    if (Do_Boundcheck()) {
      Kr_s = sol.MakeNewIteratesVectorCopy();
    }

    sol.Axpy(1.0, *IpData().trial());

    if (Do_Boundcheck()) {
      DBG_PRINT((dbg_verbosity, "Entering boundcheck"));
      // initialize
      Index new_du_size =0;
      Number* new_du_values;
      std::vector<Index> x_bound_violations_idx;
      std::vector<Number> x_bound_violations_du;
      std::vector<Index> delta_u_sort;
      bool bounds_violated;
      SmartPtr<DenseVectorSpace> delta_u_space = new DenseVectorSpace(0);
      SmartPtr<DenseVector> old_delta_u = new DenseVector(GetRawPtr(delta_u_space));
      SmartPtr<DenseVector> new_delta_u;

      bounds_violated = BoundCheck(sol, x_bound_violations_idx, x_bound_violations_du);
      while (bounds_violated) {
	Driver()->data_A()->Print(Jnlst(),J_VECTOR,J_USER1,"data_A_init");
	Driver()->data_B()->Print(Jnlst(),J_VECTOR,J_USER1,"data_B_init");
	// write new schurdata A
	dynamic_cast<IndexSchurData*>(GetRawPtr(Driver()->data_A_nonconst()))->AddData_List(x_bound_violations_idx, delta_u_sort, new_du_size, 1);
	// write new schurdata B
	dynamic_cast<IndexSchurData*>(GetRawPtr(Driver()->data_B_nonconst()))->AddData_List(x_bound_violations_idx, delta_u_sort, new_du_size, 1);
	Driver()->data_A()->Print(Jnlst(),J_VECTOR,J_USER1,"data_A");
	Driver()->data_B()->Print(Jnlst(),J_VECTOR,J_USER1,"data_B");
	Driver()->SchurBuild();
	Driver()->SchurFactorize();

	old_delta_u->Print(Jnlst(),J_VECTOR,J_USER1,"old_delta_u");
	delta_u_space = NULL; // delete old delta_u space
	delta_u_space = new DenseVectorSpace(new_du_size); // create new delta_u space
	new_delta_u = new DenseVector(GetRawPtr(ConstPtr(delta_u_space)));
	new_du_values = new_delta_u->Values();
	IpBlasDcopy(old_delta_u->Dim(), old_delta_u->Values(), 1, new_du_values, 1);
	for (UINT i=0; i<x_bound_violations_idx.size(); ++i) {
	  //	  printf("i=%d, delta_u_sort[i]=%d, x_bound_viol_du[i]=%f\n", i, delta_u_sort[i], x_bound_violations_du[i]);
	  new_du_values[delta_u_sort[i]] = x_bound_violations_du[i];
	}
	SmartPtr<IteratesVector> new_sol = sol.MakeNewIteratesVector();
	new_delta_u->Print(Jnlst(),J_VECTOR,J_USER1,"new_delta_u");

	// solve with new data_B and delta_u
	retval = Driver()->SchurSolve(&sol, ConstPtr(delta_u_long), dynamic_cast<Vector*>(GetRawPtr(new_delta_u)), Kr_s);

	sol.Axpy(1.0, *IpData().trial());

	x_bound_violations_idx.clear();
	x_bound_violations_du.clear();
	delta_u_sort.clear();
	bounds_violated = BoundCheck(sol, x_bound_violations_idx, x_bound_violations_du);
	// copy new vector in old vector ->has to be done becpause otherwise only pointers will be copied and then it makes no sense
	old_delta_u = new_delta_u->MakeNewDenseVector();
	old_delta_u->Copy(*new_delta_u);
      }
    }

    return retval;
  }
  bool StdStepCalculator::BoundCheck(IteratesVector& sol,
				     std::vector<Index>& x_bound_violations_idx,
				     std::vector<Number>& x_bound_violations_du)
  {
    DBG_START_METH("StdStepCalculator::BoundCheck", dbg_verbosity);
    DBG_ASSERT(x_bound_violations_idx.empty());
    DBG_ASSERT(x_bound_violations_du.empty());

    // find bound violations in x vector
    const Number* x_val = dynamic_cast<const DenseVector*>(GetRawPtr(IpData().curr()->x()))->Values();
    const Number* sol_val = dynamic_cast<const DenseVector*>(GetRawPtr(sol.x()))->Values();

    SmartPtr<Vector> x_L_exp = IpData().curr()->x()->MakeNew();
    SmartPtr<Vector> x_U_exp = IpData().curr()->x()->MakeNew();

    SmartPtr<Vector> x_L_comp = IpNLP().x_L()->MakeNew();
    SmartPtr<Vector> x_U_comp = IpNLP().x_U()->MakeNew();

    IpNLP().Px_L()->TransMultVector(1.0, *sol.x(), 0.0, *x_L_comp);
    IpNLP().Px_U()->TransMultVector(1.0, *sol.x(), 0.0, *x_U_comp);

    x_L_comp->Print(Jnlst(),J_VECTOR,J_USER1,"x_L_comp");
    x_U_comp->Print(Jnlst(),J_VECTOR,J_USER1,"x_U_comp");
    //    return false;

    Number* x_L_val = dynamic_cast<DenseVector*>(GetRawPtr(x_L_comp))->Values();
    Number* x_U_val = dynamic_cast<DenseVector*>(GetRawPtr(x_U_comp))->Values();

    const Number* x_L_bound = dynamic_cast<const DenseVector*>(GetRawPtr(IpNLP().x_L()))->Values();
    const Number* x_U_bound = dynamic_cast<const DenseVector*>(GetRawPtr(IpNLP().x_U()))->Values();

    for (Index i=0; i<x_L_comp->Dim(); ++i) {
      x_L_val[i] -= x_L_bound[i];
    }

    for (Index i=0; i<x_U_comp->Dim(); ++i) {
      x_U_val[i] -= x_U_bound[i];
    }

    // project back
    IpNLP().Px_L()->MultVector(1.0, *x_L_comp, 0.0, *x_L_exp);
    IpNLP().Px_U()->MultVector(1.0, *x_U_comp, 0.0, *x_U_exp);

    const Number* x_L_exp_val = dynamic_cast<DenseVector*>(GetRawPtr(x_L_exp))->Values();
    const Number* x_U_exp_val = dynamic_cast<DenseVector*>(GetRawPtr(x_U_exp))->Values();

    for (Index i=0; i<x_L_exp->Dim(); ++i) {
      if (x_L_exp_val[i]<-bound_eps_) {
	x_bound_violations_idx.push_back(i);
	x_bound_violations_du.push_back(-x_L_exp_val[i]+sol_val[i]-x_val[i]); // this is just an awkward way to compute x_bound[i] - x_curr_val[i].
      } else if (-x_U_exp_val[i]<-bound_eps_) {
	x_bound_violations_idx.push_back(i);
	x_bound_violations_du.push_back(-x_U_exp_val[i]+sol_val[i]-x_val[i]);
      }
    }

    // z_L and z_U bound violations -> These are much easier since there is no projecting back and forth
    SmartPtr<const DenseVector> z_L = dynamic_cast<const DenseVector*>(GetRawPtr(sol.z_L()));
    SmartPtr<const DenseVector> z_U = dynamic_cast<const DenseVector*>(GetRawPtr(sol.z_U()));
    z_L->Print(Jnlst(),J_VECTOR,J_USER1,"z_L_boundcheck");
    z_U->Print(Jnlst(),J_VECTOR,J_USER1,"z_U_boundcheck");
    const Number* z_L_val = z_L->Values();
    const Number* z_U_val = z_U->Values();

    SmartPtr<const DenseVector> z_L_trial = dynamic_cast<const DenseVector*>(GetRawPtr(IpData().trial()->z_L()));
    SmartPtr<const DenseVector> z_U_trial = dynamic_cast<const DenseVector*>(GetRawPtr(IpData().trial()->z_U()));
    const Number* z_L_trial_val = z_L_trial->Values();
    const Number* z_U_trial_val = z_U_trial->Values();

    // find absolute index of z_L and z_U in IteratesVector
    Index z_L_ItVec_idx = 0;
    for (Index i=0; i<4; ++i) {
      z_L_ItVec_idx += (sol.GetComp(i))->Dim();
    }
    Index z_U_ItVec_idx = z_L_ItVec_idx + sol.z_L()->Dim();

    for (Index i=0; i<z_L->Dim(); ++i) {
      if (z_L_val[i]<-bound_eps_) {
	x_bound_violations_idx.push_back(i+z_L_ItVec_idx);
	x_bound_violations_du.push_back(-z_L_trial_val[i]);
	//printf("Lower Bound Mult. no. i=%d invalid: delta_u=%f\n", i+z_L_ItVec_idx, z_L_val[i]);
      }
    }

    for (Index i=0; i<z_U->Dim(); ++i) {
      if (z_U_val[i]<-bound_eps_) {
	x_bound_violations_idx.push_back(i+z_U_ItVec_idx);
	x_bound_violations_du.push_back(-z_U_trial_val[i]);
	//printf("Upper Bound Mult. no. i=%d invalid: delta_u=%f\n", i+z_U_ItVec_idx, z_U_val[i]);
      }
    }

    //    if (x_bound_violations_idx.empty() || z_L_bound_violations_idx.empty() || z_U_bound_violations_idx.empty()) {
    if (x_bound_violations_idx.empty()) {
      return false;
    }
    else {
      return true;
    }
  }
Exemple #7
0
Number ProbingMuOracle::CalculateAffineMu(
   Number                alpha_primal,
   Number                alpha_dual,
   const IteratesVector& step
   )
{
   // Get the current values of the slack variables and bound multipliers
   SmartPtr<const Vector> slack_x_L = IpCq().curr_slack_x_L();
   SmartPtr<const Vector> slack_x_U = IpCq().curr_slack_x_U();
   SmartPtr<const Vector> slack_s_L = IpCq().curr_slack_s_L();
   SmartPtr<const Vector> slack_s_U = IpCq().curr_slack_s_U();

   SmartPtr<const Vector> z_L = IpData().curr()->z_L();
   SmartPtr<const Vector> z_U = IpData().curr()->z_U();
   SmartPtr<const Vector> v_L = IpData().curr()->v_L();
   SmartPtr<const Vector> v_U = IpData().curr()->v_U();

   SmartPtr<Vector> tmp_slack;
   SmartPtr<Vector> tmp_mult;
   SmartPtr<const Matrix> P;
   Index ncomp = 0;
   Number sum = 0.;

   // For each combination of slack and multiplier, compute the new
   // values and their dot products.

   // slack_x_L
   if( slack_x_L->Dim() > 0 )
   {
      ncomp += slack_x_L->Dim();

      P = IpNLP().Px_L();
      tmp_slack = slack_x_L->MakeNew();
      tmp_slack->Copy(*slack_x_L);
      P->TransMultVector(alpha_primal, *step.x(), 1.0, *tmp_slack);

      tmp_mult = z_L->MakeNew();
      tmp_mult->Copy(*z_L);
      tmp_mult->Axpy(alpha_dual, *step.z_L());

      sum += tmp_slack->Dot(*tmp_mult);
   }

   // slack_x_U
   if( slack_x_U->Dim() > 0 )
   {
      ncomp += slack_x_U->Dim();

      P = IpNLP().Px_U();
      tmp_slack = slack_x_U->MakeNew();
      tmp_slack->Copy(*slack_x_U);
      P->TransMultVector(-alpha_primal, *step.x(), 1.0, *tmp_slack);

      tmp_mult = z_U->MakeNew();
      tmp_mult->Copy(*z_U);
      tmp_mult->Axpy(alpha_dual, *step.z_U());

      sum += tmp_slack->Dot(*tmp_mult);
   }

   // slack_s_L
   if( slack_s_L->Dim() > 0 )
   {
      ncomp += slack_s_L->Dim();

      P = IpNLP().Pd_L();
      tmp_slack = slack_s_L->MakeNew();
      tmp_slack->Copy(*slack_s_L);
      P->TransMultVector(alpha_primal, *step.s(), 1.0, *tmp_slack);

      tmp_mult = v_L->MakeNew();
      tmp_mult->Copy(*v_L);
      tmp_mult->Axpy(alpha_dual, *step.v_L());

      sum += tmp_slack->Dot(*tmp_mult);
   }

   // slack_s_U
   if( slack_s_U->Dim() > 0 )
   {
      ncomp += slack_s_U->Dim();

      P = IpNLP().Pd_U();
      tmp_slack = slack_s_U->MakeNew();
      tmp_slack->Copy(*slack_s_U);
      P->TransMultVector(-alpha_primal, *step.s(), 1.0, *tmp_slack);

      tmp_mult = v_U->MakeNew();
      tmp_mult->Copy(*v_U);
      tmp_mult->Axpy(alpha_dual, *step.v_U());

      sum += tmp_slack->Dot(*tmp_mult);
   }

   DBG_ASSERT(ncomp > 0);

   return sum / ((Number) ncomp);
}