46 #include "AbstractLinAlgPack/src/AbstractLinAlgPack_GenPermMatrixSlice.hpp"
49 #include "Midynamic_cast_verbose.h"
51 namespace LinAlgOpPack {
55 namespace ConstrainedOptPack {
70 ,i_x_fixed_t* i_x_fixed
71 ,bnd_fixed_t* bnd_fixed
72 ,j_f_decomp_t* j_f_decomp
85 const MatrixHessianSuperBasic
86 *G_super_ptr =
dynamic_cast<const MatrixHessianSuperBasic*
>(&G);
88 if( G_super_ptr == NULL ) {
90 g,G,etaL,dL,dU,F,trans_F,f,d,nu,n_R,i_x_free,i_x_fixed,bnd_fixed
91 ,j_f_decomp,b_X,Ko,fo);
95 const MatrixHessianSuperBasic
96 &G_super = *G_super_ptr;
99 const GenPermMatrixSlice
100 &Q_R = G_super.Q_R(),
101 &Q_X = G_super.Q_X();
113 i_x_free->resize( Q_R.is_identity() ? 0: nd_R );
114 if( nd_R && !Q_R.is_identity() ) {
115 GenPermMatrixSlice::const_iterator
117 for( ; Q_itr != Q_R.end(); ++Q_itr ) {
122 (*i_x_free)[k-1] = i;
126 i_x_fixed->resize(nd_X+1);
129 GenPermMatrixSlice::const_iterator
131 for( ; Q_itr != Q_X.end(); ++Q_itr ) {
136 (*i_x_fixed)[k-1] = i;
139 (*i_x_fixed)[nd_X] = nd+1;
141 bnd_fixed->resize(nd_X+1);
144 typedef MatrixHessianSuperBasic MHSB;
145 const MHSB::bnd_fixed_t
146 &bnd_fixed_from = G_super.bnd_fixed();
148 MHSB::bnd_fixed_t::const_iterator
149 bnd_from_itr = bnd_fixed_from.begin();
150 bnd_fixed_t::iterator
151 bnd_to_itr = bnd_fixed->begin();
152 std::copy( bnd_from_itr, bnd_fixed_from.end(), bnd_to_itr );
154 (*bnd_fixed)[nd_X] =
LOWER;
156 j_f_decomp->resize(0);
164 bnd_fixed_t::const_iterator
165 bnd_itr =
const_cast<const bnd_fixed_t&
>(*bnd_fixed).begin(),
166 bnd_itr_end =
const_cast<const bnd_fixed_t&
>(*bnd_fixed).begin() + nd_X;
167 i_x_fixed_t::const_iterator
168 i_x_itr =
const_cast<const i_x_fixed_t&
>(*i_x_fixed).begin();
170 b_X_itr = b_X->begin();
171 const SpVectorSlice::element_type
173 for( ; bnd_itr != bnd_itr_end; ++bnd_itr, ++i_x_itr, ++b_X_itr ) {
178 *b_X_itr = (ele = dL.lookup_element(i))->value();
181 *b_X_itr = (ele = dU.lookup_element(i))->value();
190 *Ko = G_super.B_RR_ptr();
193 if( nd_X && G_super.B_RX_ptr().get() )
194 Vp_StMtV( &(*fo)(), -1.0, *G_super.B_RX_ptr(), G_super.B_RX_trans(), (*b_X)(1,nd_X) );
AbstractLinAlgPack::size_type size_type
int resize(OrdinalType length_in)
std::vector< EBounds > bnd_fixed_t
void V_StMtV(VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const V &V_rhs2)
v_lhs = alpha * op(M_rhs1) * V_rhs2.
QPSchurInitKKTSystemHessianFull init_kkt_full_
void V_mV(VectorMutable *v_lhs, const V &V_rhs)
v_lhs = - V_rhs.
void copy(const f_int &N, const f_dbl_prec *X, const f_int &INCX, f_dbl_prec *Y, const f_int &INCY)
T_To & dyn_cast(T_From &from)
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)
void initialize_kkt_system(const Vector &g, const MatrixOp &G, value_type etaL, const Vector *dL, const Vector *dU, const MatrixOp *F, BLAS_Cpp::Transp trans_F, const Vector *f, const Vector *d, const Vector *nu, size_type *n_R, i_x_free_t *i_x_free, i_x_fixed_t *i_x_fixed, bnd_fixed_t *bnd_fixed, j_f_decomp_t *j_f_decomp, DVector *b_X, Ko_ptr_t *Ko, DVector *fo) const
Initialize the KKT system where all variables (except the relaxation variable) are initially free and...
SparseVectorSlice< SparseElement< index_type, value_type > > SpVectorSlice
DenseLinAlgPack::VectorSliceTmpl< value_type > DVectorSlice
void initialize_kkt_system(const DVectorSlice &g, const MatrixOp &G, value_type etaL, const SpVectorSlice &dL, const SpVectorSlice &dU, const MatrixOp *F, BLAS_Cpp::Transp trans_F, const DVectorSlice *f, const DVectorSlice &d, const SpVectorSlice &nu, size_type *n_R, i_x_free_t *i_x_free, i_x_fixed_t *i_x_fixed, bnd_fixed_t *bnd_fixed, j_f_decomp_t *j_f_decomp, DVector *b_X, Ko_ptr_t *Ko, DVector *fo) const
Initialize the KKT system where the variables are initiallly fixed and free and no constraints are in...
AbstractLinAlgPack::value_type value_type
std::vector< size_type > i_x_fixed_t
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)