MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ConstrainedOptPack_QPSolverRelaxedQPOPTSOL.cpp
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 #include <assert.h>
43 
44 #include <vector>
45 
56 
57 // /////////////////////////////////////////////////////////////////
58 //
59 // This subclass uses a relaxation of the equality and inequality
60 // constraints. The mapping to the arguments of QPOPT or QPSOL
61 // is done as follows.
62 //
63 // QP formulation:
64 // ---------------
65 //
66 // min g'*d + 1/2*d'*G*d + (eta + 1/2*eta^2)*M
67 // d <: R^n
68 //
69 // s.t.
70 // etaL <= eta
71 // dL <= d <= dU
72 // eL <= op(E)*d - b*eta <= eU
73 // op(F)*d + (1 - eta) * f = 0
74 //
75 // Rearranged to :
76 // ---------------
77 //
78 // min [ g', M ] * [ d ] + 1/2 * [ d', eta ] * [ G 0 ] * [ d ]
79 // [ eta ] [ 0 M ] [ eta ]
80 //
81 // s.t. [ bL ] [ I , 0 ] [ dU ]
82 // [ etaL ] <= [ 0 , 1 ] * [ d ] <= [ inf ]
83 // [ eL ] [ op(E), -b ] [ eta ] [ eU ]
84 // [ -f ] [ op(F), -f ] [ -f ]
85 //
86 // Which maps to the QPSOL interface which is:
87 // -------------------------------------------
88 //
89 // min CVEC' * X + 1/2 * X'* H * X
90 //
91 // s.t. BL <= [ X ] <= BU
92 // [ A * X ]
93 //
94 // Which gives us:
95 //
96 // X = [ d; eta ]
97 // CVEC = [ g; M ]
98 // H = [ G, 0; 0, M ]
99 // BL = [ bL, etaL, eL, -f ]
100 // BU = [ bU, inf, eU, -f ]
101 // A = [ op(E), -b; op(F), -f ]
102 //
103 
104 namespace LinAlgOpPack {
108 }
109 
110 // ///////////////////////////////////////
111 // Members for QPSolverRelaxedQPOPTSOL
112 
113 namespace ConstrainedOptPack {
114 
116  :N_(0)
117  ,bigM_(1e+10)
118  ,use_as_bigM_(1e+10)
119  ,G_(NULL)
120 {}
121 
123 {
124  this->release_memory();
125 }
126 
127 const MatrixOp* QPSolverRelaxedQPOPTSOL::G() const
128 {
129  return G_;
130 }
131 
133 {
134  return use_as_bigM_;
135 }
136 
137 // Overridden from QPSolverRelaxed
138 
141 {
142  return qp_stats_;
143 }
144 
146 {
147  // Todo: resize to zero all the matrices and vectors
148 }
149 
152  std::ostream* out, EOutputLevel olevel, ERunTests test_what
153  ,const Vector& g, const MatrixSymOp& G
154  ,value_type etaL
155  ,const Vector* dL, const Vector* dU
156  ,const MatrixOp* E, BLAS_Cpp::Transp trans_E, const Vector* b
157  ,const Vector* eL, const Vector* eU
158  ,const MatrixOp* F, BLAS_Cpp::Transp trans_F, const Vector* f
159  ,value_type* obj_d
160  ,value_type* eta, VectorMutable* d
161  ,VectorMutable* nu
162  ,VectorMutable* mu, VectorMutable* Ed
163  ,VectorMutable* lambda, VectorMutable* Fd
164  )
165 {
166 
169 
170 #ifdef PROFILE_HACK_ENABLED
171  ProfileHackPack::ProfileTiming profile_timing( "QPSolverRelaxedQPOPTSOL::imp_solve_qp(...)" );
172 #endif
173 
174  const size_type n = d->dim();
175  const value_type inf = this->infinite_bound();
176 
177  //
178  // Map to the input arguments for QPOPT or QPSOL
179  //
180 
181  // N
182  N_ = n + 1; // With relaxation
183 
184  // NCLIN
185  n_inequ_bnds_ = ( E ? AbstractLinAlgPack::num_bounded(*eL,*eU,inf) : 0 );
186  NCLIN_ = n_inequ_bnds_ + (F ? f->dim() : 0);
187 
188  // A, BL, BU
189  A_.resize(NCLIN_,N_);
190  BL_.resize(N_+NCLIN_);
191  BU_.resize(N_+NCLIN_);
192  if(dL) {
193  VectorDenseEncap dL_de(*dL);
194  BL_(1,n) = dL_de();
195  }
196  else {
197  BL_(1,n) = -inf;
198  }
199  if(dU) {
200  VectorDenseEncap dU_de(*dU);
201  BU_(1,n) = dU_de();
202  }
203  else {
204  BU_(1,n) = -inf;
205  }
206  BL_(N_) = etaL;
207  BU_(N_) = +inf;
209  E!=NULL, std::logic_error
210  ,"Error, the QPOPT/QPSOL wrapper has not been updated for general inequalities yet!"
211  );
212 /* ToDo: Update this when needed!
213  if( E ) {
214  i_inequ_bnds_.resize(n_inequ_bnds_);
215  if( n_inequ_bnds_ < b->dim() ) {
216  // Initialize BL, BU, and A for sparse bounds on general inequalities
217  //
218  // read iterators
219  AbstractLinAlgPack::sparse_bounds_itr
220  eLU_itr( eL->begin(), eL->end(), eL->offset()
221  , eU->begin(), eU->end(), eU->offset(), inf );
222  // written iterators
223  DVector::iterator
224  BL_itr = BL_.begin() + N_,
225  BU_itr = BU_.begin() + N_;
226  ibnds_t::iterator
227  ibnds_itr = i_inequ_bnds_.begin();
228  // loop
229  for(size_type i = 1; i <= n_inequ_bnds_; ++i, ++eLU_itr, ++ibnds_itr ) {
230  TEUCHOS_TEST_FOR_EXCEPT( !( !eLU_itr.at_end() ) );
231  const size_type k = eLU_itr.indice();
232  *BL_itr++ = eLU_itr.lbound();
233  *BU_itr++ = eLU_itr.ubound();
234  *ibnds_itr = k; // Only for my record
235  // Add the corresponding row of [ op(E), -b ] to A
236  // y == A.row(i)
237  // y(1,n) = op(E')*e_k
238  DVectorSlice y = A_.row(i);
239  AbstractLinAlgPack::EtaVector e_k(k,eL->dim());
240  LinAlgOpPack::V_MtV( &y(1,n), *E, BLAS_Cpp::trans_not(trans_E), e_k() ); // op(E')*e_k
241  // y(n+1) = -b(k)
242  y(n+1) = -(*b)(k);
243  }
244  }
245  else {
246  // Initialize BL, BU and A for dense bounds on general inequalities
247  //
248  // Initialize BL(N+1:N+n_inequ_bnds), BU(N+1:N+n_inequ_bnds)
249  // and i_inequ_bnds_ = identity (only for my record, not used by QPKWIK)
250  AbstractLinAlgPack::sparse_bounds_itr
251  eLU_itr( eL->begin(), eL->end(), eL->offset()
252  , eU->begin(), eU->end(), eU->offset(), inf );
253  DVector::iterator
254  BL_itr = BL_.begin() + N_,
255  BU_itr = BU_.begin() + N_;
256  ibnds_t::iterator
257  ibnds_itr = i_inequ_bnds_.begin();
258  for(size_type i = 1; i <= n_inequ_bnds_; ++i, ++eLU_itr, ++ibnds_itr ) {
259  TEUCHOS_TEST_FOR_EXCEPT( !( !eLU_itr.at_end() ) );
260  const size_type k = eLU_itr.indice();
261  *BL_itr++ = eLU_itr.lbound();
262  *BU_itr++ = eLU_itr.ubound();
263  *ibnds_itr = k; // Only for my record
264  }
265  // A(1:n_inequ_bnds,1:n) = op(E)
266  LinAlgOpPack::assign( &A_(1,n_inequ_bnds_,1,n), *E, trans_E );
267  // A(1:n_inequ_bnds,n+1) = -b
268  LinAlgOpPack::V_StV( &A_.col(n+1)(1,n_inequ_bnds_), -1.0, *b );
269  }
270  }
271 */
273  F!=NULL, std::logic_error
274  ,"Error, the QPOPT/QPSOL wrapper has not been updated for general equalities yet!"
275  );
276 /* ToDo: Update this when needed!
277  if( F ) {
278  // BL(N+n_inequ_bnds+1:N+NCLIN) = -f
279  LinAlgOpPack::V_StV( &BL_(N_+n_inequ_bnds_+1,N_+NCLIN_), -1.0, *f );
280  // BU(N+n_inequ_bnds+1:N+NCLIN) = -f
281  LinAlgOpPack::V_StV( &BU_(N_+n_inequ_bnds_+1,N_+NCLIN_), -1.0, *f );
282  // A(n_inequ_bnds+1:NCLIN,1:n) = op(F)
283  LinAlgOpPack::assign( &A_(n_inequ_bnds_+1,NCLIN_,1,n), *F, trans_F );
284  // A(n_inequ_bnds+1:NCLIN,n+1) = -f
285  LinAlgOpPack::V_StV( &A_.col(n+1)(n_inequ_bnds_+1,NCLIN_), -1.0, *f );
286  }
287 */
288 
289  // CVEC
290  CVEC_.resize(N_);
291  CVEC_(1,n) = VectorDenseEncap(g)();
292  CVEC_(n+1) = bigM_;
293 
294  // HESS
295  G_ = &G; // That's all there is to it!
296 
297  // ISTATE
298  ISTATE_.resize(N_+NCLIN_);
299  std::fill( ISTATE_.begin(), ISTATE_.end(), 0 ); // cold start!
300  ISTATE_[n] = 1; // Make eta >= etaL active
301 
302  // X
303  X_.resize(N_);
304  X_(1,n) = VectorDenseEncap(*d)();
305  X_(n+1) = *eta;
306 
307  // AX
308  // will be resized by QPOPT but not QPSOL
309 
310  // CLAMBDA
312 
313  // LIWORK, IWORK
314  LIWORK_ = liwork(N_,NCLIN_);
315  if(static_cast<f_int>(IWORK_.size()) < LIWORK_) IWORK_.resize(LIWORK_);
316 
317  // LWORK, WORK
318  LWORK_ = lrwork(N_,NCLIN_);
319  if(static_cast<f_int>(WORK_.size()) < LWORK_) WORK_.resize(LWORK_);
320 
321  // We need to initialize some warm start information if
322  // it was given by the user!
323  bool warm_start = false;
324  if( (nu && nu->nz()) || (mu && mu->nz() ) ) {
325  // Let's try a warm start
326  if(nu) {
327  VectorDenseEncap nu_de(*nu);
328  for(int i = 1; i <= n; ++i ) {
329  if( nu_de()(i) < 0.0 )
330  ISTATE_[i-1] = 1; // Lower bound is active
331  else if( nu_de()(i) > 0.0 )
332  ISTATE_[i-1] = 2; // Upper bound is active
333  }
334  }
336  mu!=NULL, std::logic_error
337  ,"Error, the QPOPT/QPSOL wrapper has not been updated for general inequalities yet!"
338  );
339 /* ToDo: Update below when needed!
340  if(mu) {
341  const SpVectorSlice::difference_type o = mu->offset();
342  for( SpVectorSlice::const_iterator itr = mu->begin(); itr != mu->end(); ++itr ) {
343  if( itr->value() < 0.0 )
344  ISTATE_[ itr->indice() + o + n ] = 1; // Lower bound is active
345  else if( itr->value() > 0.0 )
346  ISTATE_[ itr->indice() + o + n ] = 2; // Upper bound is active
347  }
348  }
349 */
350  warm_start = true;
351  }
352 
353  //
354  // Solve the QP using QPOPT or QPSOL
355  //
356 
357  const EInform inform_return = call_qp_solver(warm_start);
358 
359  //
360  // Map from the output from QPOPT or QPSOL
361  //
362 
363  // d
364  {
365  VectorDenseMutableEncap d_de(*d);
366  d_de() = X_(1,n);
367  }
368 
369  // eta
370  *eta = X_(n+1);
371 
372  // obj_d
373  if(obj_d)
374  *obj_d = OBJ_ - (*eta) * bigM_ - 0.5 * (*eta)*(*eta) * use_as_bigM_;
375 
376  // nu
377  if(nu) {
378  VectorDenseMutableEncap nu_de(*nu);
379  nu_de() = 0.0;
380  ISTATE_t::const_iterator
381  istate_itr = ISTATE_.begin();
382  DVector::const_iterator
383  clamda_itr = CLAMDA_.begin();
384  for( size_type i = 1; i <= n; ++i, ++istate_itr, ++clamda_itr ) {
385  const f_int state = *istate_itr;
386  switch(state) {
387  case -2: // The lower bound is violated by more than feas_tol
388  case -1: // The upper bound is violated by more than feas_tol
389  // What do we do?
390  break;
391  case 0: // Within bounds by more than feas_tol
392  break;
393  case 1: // lower bound is active
394  case 2: // upper bound is active
395  case 3: // the bounds are equal and are satisfied
396  nu_de()(i) = -(*clamda_itr); // Different sign convention
397  break;
398  case 4: // Temporaraly fixed at current value
399  // What do we do?
400  break;
401  default:
402  TEUCHOS_TEST_FOR_EXCEPT(true); // Should not get here!
403  }
404  }
405  }
406 
407  // mu
409  n_inequ_bnds_!=0, std::logic_error
410  ,"Error, the QPOPT/QPSOL wrapper has not been updated for general inequalities yet!"
411  );
412 /* ToDo: Update below when needed!
413  if( n_inequ_bnds_ ) {
414  mu->resize(b->dim(),n_inequ_bnds_);
415  typedef SpVector::element_type ele_t;
416  ISTATE_t::const_iterator
417  istate_itr = ISTATE_.begin() + N_;
418  DVector::const_iterator
419  clamda_itr = CLAMDA_.begin() + N_;
420  ibnds_t::const_iterator
421  bnd_itr = i_inequ_bnds_.begin();
422  for( size_type k = 1; k <= n_inequ_bnds_; ++k, ++istate_itr, ++clamda_itr, ++bnd_itr )
423  {
424  const f_int state = *istate_itr;
425  const size_type j = *bnd_itr;
426  switch(state) {
427  case -2: // The lower bound is violated by more than feas_tol
428  case -1: // The upper bound is violated by more than feas_tol
429  // What do we do?
430  break;
431  case 0: // Within bounds by more than feas_tol
432  break;
433  case 1: // lower bound is active
434  case 2: // upper bound is active
435  case 3: // the bounds are equal and are satisfied
436  mu->add_element(ele_t(j,-(*clamda_itr))); // Different sign!
437  break;
438  case 4: // Temporaraly fixed at current value
439  // What do we do?
440  break;
441  default:
442  TEUCHOS_TEST_FOR_EXCEPT(true); // Should not get here!
443  }
444  }
445  mu->assume_sorted(true);
446  }
447  else if(E) {
448  mu->resize( eL->dim(), 0 );
449  }
450 */
451 
453  F!=NULL, std::logic_error
454  ,"Error, the QPOPT/QPSOL wrapper has not been updated for general equalities yet!"
455  );
456 /* ToDo: Update this when needed!
457 
458  // lambda
459  if( F ) {
460  LinAlgOpPack::V_StV( lambda, -1.0, CLAMDA_(N_+n_inequ_bnds_+1,N_+NCLIN_) );
461  // Validate istate
462  ISTATE_t::const_iterator
463  istate_itr = ISTATE_.begin() + N_ + n_inequ_bnds_;
464  for( size_type k = 1; k <= f->dim(); ++k, ++istate_itr ) {
465  TEUCHOS_TEST_FOR_EXCEPT( !( *istate_itr == 3 ) );
466  }
467  }
468 
469  // Ed, Fd
470  if( E && AX_.size() && eL->dim() == n_inequ_bnds_ ) {
471  if( Ed ) { // Ed = AX + b*eta
472  *Ed = AX_(1,n_inequ_bnds_);
473  if( *eta > 0.0 )
474  LinAlgOpPack::Vp_StV( Ed, *eta, *b );
475  }
476  if( Fd ) { // Fd = AX + f*eta
477  *Fd = AX_(n_inequ_bnds_+1,NCLIN_);
478  if( *eta > 0.0 )
479  LinAlgOpPack::Vp_StV( Fd, *eta, *f );
480  }
481  }
482  else {
483  if(Ed)
484  LinAlgOpPack::V_MtV( Ed, *E, trans_E, *d );
485  if(Fd)
486  LinAlgOpPack::V_MtV( Fd, *F, trans_F, *d );
487  }
488 
489 */
490 
491  //
492  // Setup the QP statistics
493  //
494 
497  switch(inform_return) {
498  case STRONG_LOCAL_MIN:
499  solution_type = QPSolverStats::OPTIMAL_SOLUTION;
500  convexity_type = QPSolverStats::CONVEX;
501  break;
502  case WEAK_LOCAL_MIN:
503  solution_type = QPSolverStats::OPTIMAL_SOLUTION;
504  convexity_type = QPSolverStats::NONCONVEX;
505  break;
506  case MAX_ITER_EXCEEDED:
507  solution_type = QPSolverStats::PRIMAL_FEASIBLE_POINT;
508  convexity_type = QPSolverStats::CONVEXITY_NOT_KNOWN;
509  break;
510  case OTHER_ERROR:
511  solution_type = QPSolverStats::SUBOPTIMAL_POINT;
512  convexity_type = QPSolverStats::CONVEXITY_NOT_KNOWN;
513  break;
514  }
516  solution_type, convexity_type
518  ,warm_start, *eta > 0.0 );
519 
520  return qp_stats_.solution_type();
521 }
522 
523 } // end namespace ConstrainedOptPack
void set_stats(ESolutionType solution_type, EConvexity convexity, int num_qp_iter, int num_adds, int num_drops, bool warm_start, bool infeasible_qp)
Initialize the statistics.
AbstractLinAlgPack::size_type size_type
virtual f_int lrwork(f_int N, f_int NCLIN) const =0
Length of real workspace.
ERunTests
Enumeration for if to run internal tests or not.
virtual const MatrixOp * G() const
Return a pointer to the matrix G to be used in the calculation of H*x by QPOPT and QPSOL...
virtual f_int liwork(f_int N, f_int NCLIN) const =0
Length of integer workspace.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void Vp_StV(VectorMutable *v_lhs, const value_type &alpha, const Vector &v_rhs)
v_lhs = alpha * v_rhs + v_lhs
EConvexity
Enumeration for the type of projected QP on output.
int resize(OrdinalType length_in)
Helper class that takes care of timing.
virtual EInform call_qp_solver(bool warm_start)=0
Solve the QP defined in the protected input data members and set the solution in the protected output...
Extract a constant DenseLinAlgPack::DVectorSlice view of a Vector object.
std::ostream * out
void Mp_StM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
void Vp_StMtV(VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta=1.0)
v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs (BLAS xGEMV)
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.
Class for storing statistics about a run of a (active set?) QP solver.
const char inf
EOutputLevel
Enumeration for the amount of output to create from solve_qp().
QPSolverStats::ESolutionType imp_solve_qp(std::ostream *out, EOutputLevel olevel, ERunTests test_what, const Vector &g, const MatrixSymOp &G, value_type etaL, const Vector *dL, const Vector *dU, const MatrixOp *E, BLAS_Cpp::Transp trans_E, const Vector *b, const Vector *eL, const Vector *eU, const MatrixOp *F, BLAS_Cpp::Transp trans_F, const Vector *f, value_type *obj_d, value_type *eta, VectorMutable *d, VectorMutable *nu, VectorMutable *mu, VectorMutable *Ed, VectorMutable *lambda, VectorMutable *Fd)
AbstractLinAlgPack::value_type value_type
virtual value_type use_as_bigM() const
Return the value of the "big M" used in the relaxation (called by QPHESS functions).
Extract a non-const DenseLinAlgPack::DVectorSlice view of a VectorMutable object. ...
Transp
TRANS.
int n
ESolutionType
Enumeration for the type of point returned from solve_qp(...).
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)