This file is indexed.

/usr/include/trilinos/ConstrainedOptPack_MatrixHessianSuperBasic.hpp is in libtrilinos-dev 10.4.0.dfsg-1ubuntu2.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
// @HEADER
// ***********************************************************************
// 
// Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
//                  Copyright (2003) Sandia Corporation
// 
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
// 
// This library is free software; you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as
// published by the Free Software Foundation; either version 2.1 of the
// License, or (at your option) any later version.
//  
// This library is distributed in the hope that it will be useful, but
// WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
// Lesser General Public License for more details.
//  
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
// USA
// Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov) 
// 
// ***********************************************************************
// @HEADER

#ifndef MATRIX_HESSIAN_SUPER_BASIC_H
#define MATRIX_HESSIAN_SUPER_BASIC_H

#include <vector>

#include "ConstrainedOptPack/src/ConstrainedOptPack_Types.hpp"
#include "AbstractLinAlgPack/src/MatrixSymWithOpFactorized.hpp"
#include "AbstractLinAlgPack/src/AbstractLinAlgPack_GenPermMatrixSlice.hpp"
#include "Miref_count_ptr.h"

namespace ConstrainedOptPack {

/** \brief Matrix class that represents a hessian matrix where only the super
 * submatrix for the super basic variables need be nonsingular.
 *
 * Given a Hessian matrix #B# and a partitioning of the variables
 * #Q = [ Q_R  Q_X ]# into free (superbasic) #Q_R'*x# and fixed (nonbasic)
 * #Q_X'*x# variables, this class represents the following matrix:
 \begin{verbatim}
  [n,n] == size(B)
  [n,n] == size(Q), Q is an othogonal permutation matrix (i.e. Q*Q' = Q'*Q = I)
  [n,n_R] == size(Q_R)
  [n,n_X] == size(Q_X)

  B = Q*Q'*B*Q*Q' = [ Q_R  Q_X ] * [ Q_R'*B*Q_R  Q_R'*B*Q_X ] * [ Q_R' ]
                                   [ Q_X'*B*Q_R  Q_X'*B*Q_X ]   [ Q_X' ]

    = [ Q_R  Q_X ] * [   B_RR     op(B_RX) ] * [ Q_R' ]
                     [ op(B_RX')   B_XX    ]   [ Q_X' ]

    = Q_R*B_RR*Q_R' + Q_R*op(B_RX)*Q_X' + Q_X*op(B_RX')*Q_R + Q_X*B_XX*Q_X'
 \end{verbatim}
 * Above, we allow the prepresentation of #op(B_RX) = Q_R'*B*Q_X# to be
 * transposed to allow for more flexibility.  Since #B_RR# and #B_XX# are
 * symmetric, we do not need to worry about transpose or not.  For this class
 * the matrix #B_RR# is required to be symmetric and nonsingular
 * (#MatrixSymWithOpFactorized# interface), but not necessarily positive definite.
 * This is the condition necessary for the Hessian when projected into the
 * active constraints at the solution for an NLP.
 * The other matrices hold no other special properties other than #B_XX# being
 * symmetric of course.
 *
 * The default compiler generated constructors are allowed and initialize the
 * matrix to uninitialized by default.
 */
class MatrixHessianSuperBasic
  : public virtual MatrixSymOp
{
public:

  /** \brief . */
  typedef Teuchos::RCP<const MatrixSymWithOpFactorized>
    B_RR_ptr_t;
  /** \brief . */
  typedef Teuchos::RCP<const MatrixOp>
    B_RX_ptr_t;
  /** \brief . */
  typedef Teuchos::RCP<const MatrixSymOp>
    B_XX_ptr_t;
  /** \brief . */
  typedef std::vector<EBounds>
    bnd_fixed_t;

  /** \brief Constructs to uninitialized.
   */
  MatrixHessianSuperBasic();

  /** \brief Initialize the matrix.
   *
   * Preconditions:\begin{itemize}
   * \item [#i_x_free != NULL#] #0 < i_x_free[l-1] <= n, l = 1...n_R# (throw ???)
   * \item [#i_x_fixed != NULL#]#0 < i_x_fixed[l-1] <= n, l = 1...n-n_R# (throw ???)
   * \item [#i_x_free != NULL && i_x_fixed != NULL#]
   *       #i_x_free[l-1] != i_x_fixed[p-1], l = 1...n_R, p = 1...n-n_X# (throw ???)
   * \item [#n_R > 0#] #B_RR_ptr.get() != NULL && B_RR_ptr->rows() == n_R && B_RR_ptr->cols() == n_R#
   *       (throw #std::invalid_argument#)
   * \item [#n_R == 0#] #B_RR_ptr.get() == NULL# (throw #std::invalid_argument#)
   * \item [#n_R < n && B_RX_ptr.get() != NULL && B_RX_trans == no_trans#]
   *       #B_RX_ptr->rows() == n_R && B_RX_ptr->cols() == n-n_R# (throw #std::invalid_argument#)
   * \item [#n_R < n && B_RX_ptr.get() != NULL && B_RX_trans == trans#]
   *       #B_RX_ptr->rows() == n-n_R && B_RX_ptr->cols() == n_R# (throw #std::invalid_argument#)
   * \item [#n_R == n#] #B_RX_ptr.get() == NULL# (throw ##std::invalid_argument#)
   * \item [#n_R < n#] #B_XX_ptr.get() != NULL && B_XX_ptr->rows() == n-n_R && B_XX_ptr->cols() == n-n_R#
   *       (throw #std::invalid_argument#)
   * \item [#n_R == n#] #B_XX_ptr.get() == NULL# (throw ##std::invalid_argument#)
   * \end{itemize}
   *
   * @param  n    [in] number of variables (used for consistency checking)
   * @param  n_R  [in] number of initially free variables (used for consistency checking)
   * @param  i_x_free
   *              [in] array (size #n_R#): #i_x_free[l-1], l = 1...n_R# defines
   *                   the matrix #Q_R# as:\\
   *                   #Q_R(:,l) = e(i_x_free[l-1]), l = 1...n_R#\\
   *                   The ordering of these indices is significant.
   *                   If #n == n_R# then #i_x_free == NULL# is allowed in which case
   *                   it is assumed to be identity.  If #n_R == 0# then of course
   *                   #i_x_free == NULL# is allowed.
   * @param  i_x_fixed
   *              [in] array (size #n_X = n - n_R#):
   *                   #i_x_fixed[l-1], l = 1...n_X# defines the matrix #Q_X# as:\\
   *                   #Q_X(:,l) = e(i_x_fixed[l-1]), l = 1...n_X#\\
   *                   The ordering of these indices is significant.
   *                   If #n_R==0# then #i_x_fixed == NULL# is allowed in which case
   *                   it is assumed to be identity.  If #n_R == n# then of course
   *                   #i_x_fixed == NULL# is allowed.
   * @param  bnd_fixed
   *              [in] array (size #n_X#):
   *                   #bnd_fixed[l-1], l = 1...n_X# defines what bound the variables
   *                   in #i_x_fixed[l-1], l = 1...n_X# are fixed at: #LOWER#, #UPPER#
   *                   or #EQUALITY#.  If #n_R == n# then of course
   *                   #bnd_fixed == NULL# is allowed.
   * @param  B_RR_ptr
   *              [in] Smart pointer to matrix #B_RR# (size #n_R x n_R#) for the
   *                   free (super basic) variables. #B_RR_ptr.get() != NULL#
   *                   must be true if #n_R > # or an exception will be thrown.
   *                   if #n_R == 0# then #B_RR_ptr.get() == NULL# may be true.
   * @param  B_RX_ptr
   *              [in] Smart pointer to matrix #B_RX# (size #n_R x n_X#
   *                   if #B_RX_trans==no_trans#  or #n_X x n_R# if #B_RX_trans==trans#)
   *                   for the  cross terms of free (super basic) and fixed (nonbasic)
   *                   variables.  It is allowed for #B_RX_ptr.get() == NULL#.
   * @param  B_RX_trans
   *              [in] Determines if op(B_RX) = B_RX (#no_trans#) or op(B_RX) = B_RX'
   *                   (#trans#).  Ignored if #n_R == n#.
   * @param  B_XX_ptr
   *              [in] Smart pointer to matrix B_XX (size #n_X x n_X#) for the
   *                   fixed (nonbasic) variables. #B_XX_ptr.get() != NULL#
   *                   must be true if #n_R < n# or an exception will be thrown.
   *                   If #n_R == n# then #B_XX_ptr.get() == NULL# may be true.
   */
  virtual void initialize(
    size_type            n
    ,size_type           n_R
    ,const size_type     i_x_free[]
    ,const size_type     i_x_fixed[]
    ,const EBounds       bnd_fixed[]
    ,const B_RR_ptr_t&   B_RR_ptr
    ,const B_RX_ptr_t&   B_RX_ptr
    ,BLAS_Cpp::Transp    B_RX_trans
    ,const B_XX_ptr_t&   B_XX_ptr
    );

  /** @name Provide access to constituent members */
  //@{

  /** \brief . */
  const GenPermMatrixSlice& Q_R() const;
  /** \brief . */
  const GenPermMatrixSlice& Q_X() const;
  /** \brief . */
  const bnd_fixed_t& bnd_fixed() const;
  /** \brief . */
  const B_RR_ptr_t& B_RR_ptr() const;
  /** \brief . */
  const B_RX_ptr_t& B_RX_ptr() const;
  /** \brief . */
  BLAS_Cpp::Transp B_RX_trans() const;
  /** \brief . */
  const B_XX_ptr_t& B_XX_ptr() const;

  //@}

  /** @name Overridden from Matrix */
  //@{

  /// 
  size_type rows() const;

  //@}

  /** @name Overridden from MatrixOp */
  //@{

  /** \brief . */
  void Vp_StMtV(DVectorSlice* vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
    , const DVectorSlice& vs_rhs2, value_type beta) const;
  /** \brief . */
  void Vp_StMtV(DVectorSlice* vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
    , const SpVectorSlice& sv_rhs2, value_type beta) const;
  /** \brief . */
  void Vp_StPtMtV(DVectorSlice* vs_lhs, value_type alpha
    , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
    , BLAS_Cpp::Transp M_rhs2_trans
    , const DVectorSlice& sv_rhs3, value_type beta) const;
  /** \brief . */
  value_type transVtMtV(const SpVectorSlice& sv_rhs1, BLAS_Cpp::Transp trans_rhs2
    , const SpVectorSlice& sv_rhs3) const ;

  //@}

protected:

  /** \brief . */
  void assert_initialized() const;

private:

  // ///////////////////////////////////
  // Private types

  typedef std::vector<size_type>		row_i_t;
  typedef std::vector<size_type>		col_j_t;
  
  // ///////////////////////////////////
  // Private data members

  size_type               n_;
  size_type               n_R_;
  GenPermMatrixSlice		Q_R_;        // Sorted by row
  row_i_t					Q_R_row_i_;
  col_j_t					Q_R_col_j_;
  GenPermMatrixSlice		Q_X_;        // Sorted by row
  row_i_t					Q_X_row_i_;
  col_j_t					Q_X_col_j_;
  bnd_fixed_t             bnd_fixed_;
  B_RR_ptr_t              B_RR_ptr_;
  B_RX_ptr_t              B_RX_ptr_;
  BLAS_Cpp::Transp        B_RX_trans_;
  B_XX_ptr_t              B_XX_ptr_;

}; // end class MatrixHessianSuperBasic

// ////////////////////////////////////////////
// Inline members for MatrixHessianSuperBasic

inline
const GenPermMatrixSlice& MatrixHessianSuperBasic::Q_R() const
{
  assert_initialized();
  return Q_R_;
}

inline
const GenPermMatrixSlice& MatrixHessianSuperBasic::Q_X() const
{
  assert_initialized();
  return Q_X_;
}

inline
const MatrixHessianSuperBasic::bnd_fixed_t&
MatrixHessianSuperBasic::bnd_fixed() const
{
  return bnd_fixed_;
}

inline
const MatrixHessianSuperBasic::B_RR_ptr_t&
MatrixHessianSuperBasic::B_RR_ptr() const
{
  assert_initialized();
  return B_RR_ptr_;
}

inline
const MatrixHessianSuperBasic::B_RX_ptr_t&
MatrixHessianSuperBasic::B_RX_ptr() const{
  assert_initialized();
  return B_RX_ptr_;
}

inline
BLAS_Cpp::Transp MatrixHessianSuperBasic::B_RX_trans() const
{
  assert_initialized();
  return B_RX_trans_;
}

inline
const MatrixHessianSuperBasic::B_XX_ptr_t&
MatrixHessianSuperBasic::B_XX_ptr() const
{
  assert_initialized();
  return B_XX_ptr_;
}

} // end namespace ConstrainedOptPack

#endif // MATRIX_HESSIAN_SUPER_BASIC_H