MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions
Collaboration diagram for Reduction operations:

Functions

value_type AbstractLinAlgPack::max_element (const Vector &v)
 Compute the maximum element in a vector. More...
 
std::pair< value_type, value_type > AbstractLinAlgPack::max_near_feas_step (const Vector &x, const Vector &d, const Vector &xl, const Vector &xu, value_type max_bnd_viol)
 Computes the maximum positive and negative step that can be taken that are within the relaxed bounds. More...
 
value_type AbstractLinAlgPack::max_rel_step (const Vector &x, const Vector &d)
 Computes the maximum relative step of x = x + d. More...
 
value_type AbstractLinAlgPack::fraction_to_boundary (const value_type tau, const Vector &x, const Vector &d, const Vector &xl, const Vector &xu)
 
value_type AbstractLinAlgPack::fraction_to_zero_boundary (const value_type tau, const Vector &x, const Vector &d)
 
size_type AbstractLinAlgPack::num_bounded (const Vector &xl, const Vector &xu, value_type inf_bound)
 Count the number of finitly bounded elements in xl <= x <= xu. More...
 
value_type AbstractLinAlgPack::log_bound_barrier (const Vector &x, const Vector &xl, const Vector &xu)
 Computes the log barrier term: More...
 
value_type AbstractLinAlgPack::combined_nu_comp_err (const Vector &v, const Vector &x, const Vector &xl, const Vector &xu)
 Computes an estimate of the complementarity error. More...
 
value_type AbstractLinAlgPack::combined_nu_comp_err_lower (const Vector &v, const Vector &x, const Vector &xl)
 Computes an estimate of the complementarity error when only the lower bounds are non-infinite. More...
 
value_type AbstractLinAlgPack::combined_nu_comp_err_upper (const Vector &v, const Vector &x, const Vector &xu)
 Computes an estimate of the complementarity error when only the upper bounds are non-infinite. More...
 
value_type AbstractLinAlgPack::IP_comp_err_with_mu (const value_type mu, const value_type inf_bound, const Vector &x, const Vector &xl, const Vector &xu, const Vector &vl, const Vector &vu)
 Computes the complementarity error for a primal/dual interior point algorithm using inf norm. More...
 
bool AbstractLinAlgPack::max_inequ_viol (const AbstractLinAlgPack::Vector &v, const AbstractLinAlgPack::Vector &vL, const AbstractLinAlgPack::Vector &vU, AbstractLinAlgPack::size_type *max_viol_i, AbstractLinAlgPack::value_type *max_viol, AbstractLinAlgPack::value_type *v_i, int *bnd_type, AbstractLinAlgPack::value_type *vLU_i)
 Compute the maximum violation from a set of inequality constraints vL <= v <= vU. More...
 

Detailed Description

Function Documentation

AbstractLinAlgPack::value_type AbstractLinAlgPack::max_element ( const Vector &  v)

Compute the maximum element in a vector.

Returns:

max{ v(i), i = 1...n }

Definition at line 173 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

std::pair< AbstractLinAlgPack::value_type, AbstractLinAlgPack::value_type > AbstractLinAlgPack::max_near_feas_step ( const Vector &  x,
const Vector &  d,
const Vector &  xl,
const Vector &  xu,
value_type  max_bnd_viol 
)

Computes the maximum positive and negative step that can be taken that are within the relaxed bounds.

This function computes and returns the maximum (in magnitude) postive (return.first >= 0.0) and negative (return.second <= 0.0) steps u that can be taken such that the relaxed bounds:

xl - max_bnd_viol <= x + u * d <= xu - max_bnd_viol

are strictly satisfied.

If return.first < 0.0 then this is a flag that x is not in the relaxed bounds to begin with. In this case return.second has no meaning.

Definition at line 184 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::max_rel_step ( const Vector &  x,
const Vector &  d 
)

Computes the maximum relative step of x = x + d.

return = max{ |d|/(1.0+|x(i)|), for i = 1...n }

Definition at line 204 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::fraction_to_boundary ( const value_type  tau,
const Vector &  x,
const Vector &  d,
const Vector &  xl,
const Vector &  xu 
)

Computes alpha_max by the fraction to boundary rule.

ToDo: Finish documentation!

Definition at line 220 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::fraction_to_zero_boundary ( const value_type  tau,
const Vector &  x,
const Vector &  d 
)

Computes alpha_max by the fraction to boundary rule assuming only a lower bound of zero

ToDo: Finish documentation!

Definition at line 240 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::size_type AbstractLinAlgPack::num_bounded ( const Vector &  xl,
const Vector &  xu,
value_type  inf_bound 
)

Count the number of finitly bounded elements in xl <= x <= xu.

ToDo: Finish documentation!

Definition at line 258 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::log_bound_barrier ( const Vector &  x,
const Vector &  xl,
const Vector &  xu 
)

Computes the log barrier term:

sum{ log( x(i) - xl(i) ) + log( xu(i) - x(i) ) , for i = 1...n }

Definition at line 275 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::combined_nu_comp_err ( const Vector &  v,
const Vector &  x,
const Vector &  xl,
const Vector &  xu 
)

Computes an estimate of the complementarity error.

  for every i...
     comp_err = max(comp_err, v(i)*(xu(i)-x(i), -v(i)*(x(i)-xl(i))));

Definition at line 294 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::combined_nu_comp_err_lower ( const Vector &  v,
const Vector &  x,
const Vector &  xl 
)

Computes an estimate of the complementarity error when only the lower bounds are non-infinite.

  for every i...
     comp_err = max(comp_err, v(i)*(xl(i)-x(i))));

NOTE: equivalent to 
     comp_err = max(comp_err, -v(i)*(x(i)-xl(i))));

Definition at line 313 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::combined_nu_comp_err_upper ( const Vector &  v,
const Vector &  x,
const Vector &  xu 
)

Computes an estimate of the complementarity error when only the upper bounds are non-infinite.

  for every i...
     comp_err = max(comp_err, v(i)*(xu(i)-x(i))));

Definition at line 332 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

AbstractLinAlgPack::value_type AbstractLinAlgPack::IP_comp_err_with_mu ( const value_type  mu,
const value_type  inf_bound,
const Vector &  x,
const Vector &  xl,
const Vector &  xu,
const Vector &  vl,
const Vector &  vu 
)

Computes the complementarity error for a primal/dual interior point algorithm using inf norm.

  for every i...
     comp_err = max(comp_err, 
                  ( fabs(vu(i)*(xu(i)-x(i)) - mu) ),
          ( fabs(vl(i)*(x(i)-xl(i)) - mu) )
          );

Definition at line 350 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.

bool AbstractLinAlgPack::max_inequ_viol ( const AbstractLinAlgPack::Vector v,
const AbstractLinAlgPack::Vector vL,
const AbstractLinAlgPack::Vector vU,
AbstractLinAlgPack::size_type max_viol_i,
AbstractLinAlgPack::value_type max_viol,
AbstractLinAlgPack::value_type v_i,
int *  bnd_type,
AbstractLinAlgPack::value_type vLU_i 
)

Compute the maximum violation from a set of inequality constraints vL <= v <= vU.

Parameters
v[in] Inequality value vector.
vL[in] Lower inequality bounds (may be -infinity (i.e. very large negative number))
vU[in] Upper inequality bounds (may be +infinity (i.e. very large positive number))
max_viol_i[out] Gives the index of the inequality with the maximum (scaled) violation. If *max_viol_i == 0 on output then no inequality was violated.
max_viol[out] The maximum (scaled violation). Only significant if *max_viol_i > 0.
v_i[out] Set to v.get_ele(*max_viol_i). Only significant if *max_viol_i > 0.
bnd_type[out] The type of bound with the maximum violation.
  • -1 : LOWER
  • 0 : EQUALITY
  • +1 : UPPER
Only significant if *max_viol_i > 0.
vLU_i[out] Set to:
  • vL.get_ele(*max_viol_i) if *bnd_type <= 0
  • vU.get_ele(*max_viol_i) if *bnd_type > 0
Only significant if *max_viol_i > 0.
Returns
Returns true if some constraint was violated.

Preconditions:

  • ToDo: Spell these out!

Postconditions:

  • ToDo: Spell these out!

In order to make the result unique if more than one inequality vL(i) <= v(i) <= vL(i) have the same maximum violation then the inequality with the lowest i is returned.

Definition at line 372 of file AbstractLinAlgPack_VectorAuxiliaryOps.cpp.