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; }