MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_VectorAuxiliaryOps.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef ABSTRACT_LINALG_PACK_VECTOR_AUXILIARY_OPS_H
43 #define ABSTRACT_LINALG_PACK_VECTOR_AUXILIARY_OPS_H
44 
45 #include <utility>
46 
48 
49 namespace AbstractLinAlgPack {
50 
54 
57 
66 value_type max_element( const Vector& v );
67 
83 std::pair<value_type,value_type>
85  const Vector& x, const Vector& d
86  ,const Vector& xl, const Vector& xu
87  ,value_type max_bnd_viol
88  );
89 
98  const Vector& x, const Vector& d
99  );
100 
102 
107  const value_type tau,
108  const Vector &x,
109  const Vector &d,
110  const Vector &xl,
111  const Vector &xu
112  );
113 
115 
121  const value_type tau,
122  const Vector &x,
123  const Vector &d
124  );
125 
131  const Vector& xl, const Vector& xu
132  ,value_type inf_bound );
133 
142  const Vector &x
143  ,const Vector &xl
144  ,const Vector &xu
145  );
146 
156  const Vector &v
157  ,const Vector &x
158  ,const Vector &xl
159  ,const Vector &xu
160  );
161 
162 
177  const Vector &v
178  ,const Vector &x
179  ,const Vector &xl
180  );
181 
192  const Vector &v
193  ,const Vector &x
194  ,const Vector &xu
195  );
196 
197 
211  const value_type mu
212  ,const value_type inf_bound
213  ,const Vector &x
214  ,const Vector &xl
215  ,const Vector &xu
216  ,const Vector &vl
217  ,const Vector &vu
218  );
219 
265 bool max_inequ_viol(
267  ,const AbstractLinAlgPack::Vector &vL
268  ,const AbstractLinAlgPack::Vector &vU
269  ,AbstractLinAlgPack::size_type *max_viol_i
272  ,int *bnd_type
274  );
275 
277 
280 
292 void force_in_bounds( const Vector& xl, const Vector& xu, VectorMutable* x );
293 
299  const value_type rel_push,
300  const value_type abs_push,
301  const Vector &xl,
302  const Vector &xu,
303  VectorMutable *x
304  );
305 
313 void inv_of_difference(
314  const value_type alpha
315  ,const Vector &v0
316  ,const Vector &v1
317  ,VectorMutable *z
318  );
319 
328  const Vector &xl
329  ,const value_type inf_bound_limit
330  ,VectorMutable *vl
331  );
332 
341  const Vector &xu
342  ,const value_type inf_bound_limit
343  ,VectorMutable *vu
344  );
345 
353  const value_type mu,
354  const Vector &invXl,
355  const Vector &vl,
356  const Vector &d_k,
357  VectorMutable *dvl
358  );
359 
368  const value_type mu,
369  const Vector &invXu,
370  const Vector &vu,
371  const Vector &d_k,
372  VectorMutable *dvu
373  );
374 
375 
384 void ele_wise_sqrt(
385  VectorMutable* z
386  );
387 
395 void max_vec_scalar(
396  value_type min_ele
397  ,VectorMutable *y
398  );
399 
407 void max_abs_vec_scalar(
408  value_type min_ele
409  ,VectorMutable *y
410  );
411 
413 
415 
416 } // end namespace AbstractLinAlgPack
417 
418 #endif // ABSTRACT_LINALG_PACK_VECTOR_AUXILIARY_OPS_H
void max_vec_scalar(value_type min_ele, VectorMutable *y)
Take the maximum value of the vector elements and a scalar.
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
value_type 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...
void inv_of_difference(const value_type alpha, const Vector &v0, const Vector &v1, VectorMutable *z)
Computes the inverse of the difference between two vectors.
void force_in_bounds_buffer(const value_type rel_push, const value_type abs_push, const Vector &xl, const Vector &xu, VectorMutable *x)
Force a vector sufficiently within bounds according to a specified absolute and relative buffer...
void upperbound_multipliers_step(const value_type mu, const Vector &invXu, const Vector &vu, const Vector &d_k, VectorMutable *dvu)
Calculates the multiplier step for the upper bounds.
value_type max_element(const Vector &v)
Compute the maximum element in a vector.
void max_abs_vec_scalar(value_type min_ele, VectorMutable *y)
Take the maximum value of the absolute vector elements and a scalar.
value_type max_rel_step(const Vector &x, const Vector &d)
Computes the maximum relative step of x = x + d.
value_type combined_nu_comp_err(const Vector &v, const Vector &x, const Vector &xl, const Vector &xu)
Computes an estimate of the complementarity error.
bool 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. ...
size_type num_bounded(const Vector &xl, const Vector &xu, value_type inf_bound)
Count the number of finitly bounded elements in xl <= x <= xu.
value_type fraction_to_boundary(const value_type tau, const Vector &x, const Vector &d, const Vector &xl, const Vector &xu)
value_type fraction_to_zero_boundary(const value_type tau, const Vector &x, const Vector &d)
void ele_wise_sqrt(VectorMutable *z)
Calculates the sqrt of each element in the vector Pre Condition: all elements of z must be positive...
value_type log_bound_barrier(const Vector &x, const Vector &xl, const Vector &xu)
Computes the log barrier term:
std::pair< value_type, value_type > 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...
void force_in_bounds(const Vector &xl, const Vector &xu, VectorMutable *x)
Force a vector in bounds.
void correct_lower_bound_multipliers(const Vector &xl, const value_type inf_bound_limit, VectorMutable *vl)
Corrects the lower bound multipliers with infinite bounds.
void lowerbound_multipliers_step(const value_type mu, const Vector &invXl, const Vector &vl, const Vector &d_k, VectorMutable *dvl)
Calculates the multiplier step for lower bounds.
void correct_upper_bound_multipliers(const Vector &xu, const value_type inf_bound_limit, VectorMutable *vu)
Corrects the upper bound multipliers with infinite bounds.
value_type 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.
value_type 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.