This file is indexed.

/usr/include/trilinos/ConstrainedOptPack_DecompositionSystem.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
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
// @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 DECOMPOSITION_SYSTEM_H
#define DECOMPOSITION_SYSTEM_H

#include <stdexcept>

#include "ConstrainedOptPack_Types.hpp"
#include "AbstractLinAlgPack_VectorSpace.hpp"

namespace ConstrainedOptPack {

/** \brief This class abstracts a decomposition choice for the quasi-range
 * space \a Y and null space \a Z matrices for a linearly independent
 * set of columns of \a Gc.
 *
 * <tt>Gc = [ Gc(:,equ_decomp),  Gc(:,equ_undecomp) ]</tt>
 *
 * where \c Gc is <tt>n x m</tt>, \c Gc(:,equ_decomp) is <tt>n x r</tt> and
 * \c Gc(:,equ_undecomp) is <tt>n x (m - r)</tt>.
 *
 * Note that the columns in <tt>Gc(:,equ_undecomp)</tt> may be
 * linearly dependent with the columns in <tt>Gc(:,equ_undecomp)</tt>
 * or they may just be undecomposed linearly independent equality
 * constraints.
 *
 * The decomposition formed by subclasses must have the properties:
 \verbatim
   Z s.t. Gc(:,equ_decomp)' * Z = 0
   Y s.t. [Z  Y] is nonsingular
   R = Gc(:,equ_decomp)' * Y is nonsingular
   Uz = Gc(:,equ_undecomp)' * Z
   Uy = Gc(:,equ_undecomp)' * Y
 \endverbatim
 *
 * The matrix factory objects returned by ??? are ment to have a lifetime that is
 * independent of \c this.
 *
 * The decomposition matrices \c Z, \c Y, \c R, \c Uz and \c Uy which
 * are updated in <tt>this->update_decomp()</tt> must be completely independent from
 * \c this and from each other and \c Gc that they based on.  For example,
 * Once \c update_decomp() is called, \c this, \c Gc can be destroyed and
 * the behaviors of the decomposition matrices must not be altered.  In this respect
 * the <tt>%DecompositionSystem</tt> interface is really nothing more than a "Strategy"
 * interface (with some state data of course) for computing range/null decompositions.
 * This gives the client great flexibility in how the decomposition matrices are used. 
 *
 * ToDo: Finish documentation!
 */
class DecompositionSystem {
public:

  /** @name Public types */
  //@{

  /** \brief . */
  typedef Teuchos::RCP<
    const Teuchos::AbstractFactory<MatrixOpNonsing> >    mat_nonsing_fcty_ptr_t;
  /** \brief . */
  typedef Teuchos::RCP<
    const Teuchos::AbstractFactory<MatrixOp> >               mat_fcty_ptr_t;
  /** \brief . */
  class SingularDecomposition : public std::logic_error
  {public: SingularDecomposition(const std::string& what_arg) : std::logic_error(what_arg) {}};
  /** \brief . */
  class InvalidMatrixType : public std::logic_error
  {public: InvalidMatrixType(const std::string& what_arg) : std::logic_error(what_arg) {}};
  /** \brief . */
  class TestFailed : public std::runtime_error
  {public: TestFailed(const std::string& what_arg) : std::runtime_error(what_arg) {}};
  /// Enumeration for the amount of output to create from <tt>update_decomp()</tt>.
  enum EOutputLevel {
    PRINT_NONE          = 0,
    PRINT_BASIC_INFO    = 1,
    PRINT_MORE_INFO     = 2,
    PRINT_VECTORS       = 3,
    PRINT_EVERY_THING   = 4
    };
  /// Enumeration for if to run internal tests or not.
  enum ERunTests { RUN_TESTS, NO_TESTS };
  /** \brief . */
  enum EMatRelations { MATRICES_INDEP_IMPS, MATRICES_ALLOW_DEP_IMPS };

  //@}

  /** \brief . */
  virtual ~DecompositionSystem() {}

  /** @name Dimensionality of the decomposition */
  //@{

  /** \brief Return the number of rows in \c Gc.
   *
   * Postconditions:<ul>
   * <li> <tt>n > m</tt>
   * </ul>
   *
   * The default implementation returns
   * <tt>this->space_range()->dim() + this->space_null()->dim()</tt>.
   */
  virtual size_type n() const;

  /** \brief Return the number of columns in \c Gc.
   *
   * Postconditions:<ul>
   * <li> <tt>m > 0</tt>
   * </ul>
   */
  virtual size_type m() const = 0;

  /** \brief Returns the rank of \c Gc(:,equ_decomp()).
   *
   * Postconditions:<ul>
   * <li> <tt>r =< m</tt>
   * </ul>
   *
   * The default implementation returns
   * <tt>this->space_range()->dim()</tt>.
   */
  virtual size_type r() const;

  /** \brief Returns the range of the decomposed equalities.
   *
   * The default implementation returns <tt>Range1D(1,this->r())</tt>.
   */
  virtual Range1D equ_decomp() const;

  /** \brief Returns the range of the undecomposed equalities.
   *
   * The default implementation returns <tt>Range1D(this->r()+1,this->m())</tt>
   * or <tt>Range1D::Invalid</tt> if <tt>this->r() == this->m()<tt>
   */
  virtual Range1D equ_undecomp() const;

  //@}

  /** @name Range and null vector spaces */
  //@{

  /** \brief Return a \c VectorSpace object for the range space.
   *
   * Postconditions:<ul>
   * <li> <tt>return.get() != NULL</tt>
   * <li> <tt>return->dim() == this->r()</tt>
   * </ul>
   */
  virtual const VectorSpace::space_ptr_t space_range() const = 0;

  /** \brief Return a \c VectorSpace object for the range space.
   *
   * Postconditions:<ul>
   * <li> <tt>return.get() != NULL</tt>
   * <li> <tt>return->dim() == this->n() - this->r()</tt>
   * </ul>
   */
  virtual const VectorSpace::space_ptr_t space_null() const = 0;

  //@}

  /** @name Matrix factories */
  //@{

  /** \brief Return a matrix factory object for <tt>Z</tt>
   */
  virtual const mat_fcty_ptr_t factory_Z() const = 0;

  /** \brief Return a matrix factory object for <tt>Y</tt>
   */
  virtual const mat_fcty_ptr_t factory_Y() const = 0;

  /** \brief Return a matrix factory object for <tt>R</tt>.
   */
  virtual const mat_nonsing_fcty_ptr_t factory_R() const = 0;
  
  /** \brief Return a matrix factory object for <tt>Uz</tt>
   */
  virtual const mat_fcty_ptr_t factory_Uz() const = 0;

  /** \brief Return a matrix factory object for <tt>Uy</tt>
   */
  virtual const mat_fcty_ptr_t factory_Uy() const = 0;

  //@}

  /** @name Update range/null decomposition */
  //@{

  /** \brief Creates the range/null decomposition for <tt>Gc(:,equ_decomp)'</tt>.
   *
   * The decomposition is based on the linearly independent columns \c Gc(:,equ_decomp)
   * of \c Gc
   *
   * <tt>Gc = [ Gc(:,equ_decomp),  Gc(:,equ_undecomp) ]</tt>
   *
   * Specifically this operation finds the matrices:
   \verbatim
   Z s.t. Gc(:,equ_deomp)' * Z = 0
   Y s.t. [Z  Y] is nonsingular
   R = Gc(:,equ_decomp)' * Y is nonsingular
   Uz = Gc(:,equ_undecomp)' * Z
   Uy = Gc(:,equ_undecomp)' * Y
   \endverbatim
   * If there is some problem creating the decomposition then exceptions
   * with the base class \c std::exception may be thrown.  The meaning
   * of these exceptions are more associated with the subclasses
   * that implement this operation.
   *
   * The concrete types for <tt>Gc</tt>, <tt>Z</tt>, <tt>Y</tt>,
   * <tt>Uz</tt> and <tt>Uy</tt> must be compatable with
   * the concrete implementation of \c this or an <tt>InvalidMatrixType</tt> exeption
   * will be thrown..
   *
   * @param  out [out] If <tt>out!=NULL</tt> then output is printed to this stream
   *             depending on the value of \c olevel.
   * @param olevel
   *             [in] Determines the amount of output to print to \c *out.
   *             The exact type of output is determined by the implementing
   *             subclass but here is the sugguested behavior:<ul>
   *             <li> \c PRINT_NONE : Don't print anything (same as <tt>out==NULL</tt>).
   *             <li> \c PRINT_BASIC_INFO : Only print basic information about
   *                  how the decomposition is formed.
   *                  Amount of output = \c O(1).
   *             <li> \c PRINT_VECTORS : Prints out important vectors computed
   *                  durring the computations (usually only durring testing).
   *                  This level is only useful for debugging.
   *                  Amount of output = \c O(n).
   *             <li> \c PRINT_EVERY_THING : Print out nearly every important
   *                  quantity that is computed (except for the output matrices
   *                  themselves, clients can do that) while the matrices are being
   *                  formed or tests are being conducted.  This level is only
   *                  useful for debugging.
   *                  Amount of output = <tt>O(m*n)</tt>
   *             </ul>
   * @param test_what
   *             [in] Determines if internal validation tests are performed.
   *             The post conditions for the output matrices are not checked
   *             internally.  This is something that client can (and should)
   *             do independently (see \c DecompositionSystemTester).  Values:<ul>
   *             <li> \c RUN_TESTS : As many validation/consistency tests
   *                  are performed internally as possible.  If a test
   *                  fails then a \c TestFailed execption will be thrown.
   *                  The subclasses determine what the tests are and
   *                  what failing a test means.
   *             <li> \c NO_TEST : No tests are performed internally.  This is
   *                  to allow the fastest possible execution.
   *             </ul>
   *             If a test fails, then a \c TestFailed exception will be thrown with
   *             a helpful error message.
   * @param  Gc  [in] The matrix for which the range/null decomposition is defined.
   * @param  Z   [out] On output represents the <tt>n x (n-r)</tt> null space	matrix such that
   *             <tt>Gc(:,equ_decomp) * Z == 0</tt>.  This matrix object must have been created
   *             by <tt>this->factory_Z()->create()</tt>.
   * @param  Y   [out] On output represents the <tt>n x r</tt> range space matrix	such that
   *             <tt>[ Y  Z ]</tt> is nonsingular.  This matrix object must have been created
   *             by <tt>this->factory_Y()->create()</tt>.
   * @param  R   [out] On output represents the nonsingular <tt>r x r</tt> matrix <tt>Gc(:,equ_decomp) * Y</tt>.
   *             This matrix object must have been created by <tt>this->factory_R()->create()</tt>.
   * @param  Uz  [in/out] If <tt>Uz != NULL</tt> (<tt>this->m() > this->r()</tt> only) then on output
   *             <tt>*Uz</tt> represents the <tt>(m-r) x (n-r)</tt> matrix <tt>Gc(:,equ_undecomp) * Z</tt>.
   *             If <tt>this->m() == this->r()</tt> then <tt>Uz == NULL</tt> must be true.
   *             If <tt>Uz!=NULL</tt>, then this matrix object must have been created by
   *             <tt>this->factory_Uz()->create()</tt>.
   * @param  Uy  [in/out] If <tt>Uy != NULL</tt> (<tt>this->m() > this->r()</tt> only) then on output
   *             <tt>*Uy</tt> represents the <tt>(m-r) x r</tt> matrix <tt>Gc(:,equ_undecomp) * Y</tt>.
   *             If <tt>this->m() == this->r()</tt> then <tt>Uy == NULL</tt> must be true.
   *             If <tt>Uy!=NULL</tt>, then this matrix object must have been created by
   *             <tt>this->factory_Uy()->create()</tt>.
   * @param mat_rel
   *             [in] Determines if the ouput matrix objects must be completely independent or not.
   *             <ul>
   *             <li> MATRICES_INDEP_IMPS: The matrix objects must have independent implementations (default).
   *             <li> MATRICES_ALLOW_DEP_IMPS: The matrix objects can have implementation dependencies.
   *             </ul>
   *
   * Preconditions:<ul>
   * <li> <tt>Gc.rows() == this->n()</tt> (throw \c std::invalid_argument)
   * <li> <tt>Gc.cols() == this->m()</tt> (throw \c std::invalid_argument)
   * <li> [<tt>this->m() == this->r()</tt>] <tt>Uz == NULL</tt> (throw \c std::invalid_argument)
   * <li> [<tt>this->m() == this->r()</tt>] <tt>Uy == NULL</tt> (throw \c std::invalid_argument)
   * <li> <tt>Z!=NULL || Y!=NULL || R!=NULL || Uz!=NULL || Uy!=NULL</tt>
   *      (throw \c std::invalid_argument)
   * </ul>
   *
   * Postconditions:<ul>
   * <li> [<tt>Z != NULL</tt>] <tt>Z.space_cols().is_compatible(Gc.space_cols()) == true)</tt>
   * <li> [<tt>Z != NULL</tt>] <tt>Z.space_rows().is_compatible(*space_null()) == true</tt>
   * <li> [<tt>Z != NULL</tt>] <tt>Gc(:,equ_decomp())' * Z == 0</tt>
   * <li> [<tt>Y != NULL</tt>] <tt>Y.space_cols().is_compatible(Gc.space_cols()) == true)</tt>
   * <li> [<tt>Y != NULL</tt>] <tt>Y.cols() == this->r()</tt>
   * <li> [<tt>Y != NULL</tt>] <tt>[ Y  Z ]</tt> is nonsingular
   * <li> [<tt>R != NULL</tt>] <tt>R->space_cols().is_compatible(*Gc.space_cols()->sub_space(equ_decomp())) == true</tt>
   * <li> [<tt>R != NULL</tt>] <tt>R->space_rows().is_compatible(*space_range()) == true</tt>
   * <li> [<tt>R != NULL</tt>] <tt>R == Gc(:,equ_decomp())'*Y</tt>
   * <li> [<tt>Uz != NULL</tt>] <tt>Uz.space_cols().is_compatible(*Gc.space_rows()->sub_space(equ_undecomp())) == true</tt>
   * <li> [<tt>Uz != NULL</tt>] <tt>Uz.space_rows().is_compatible(*space_null()) == true</tt>
   * <li> [<tt>Uz != NULL</tt>] <tt>Uz == Gc(:,equ_undecomp())'*Z</tt>
   * <li> [<tt>Uy != NULL</tt>] <tt>Uy.space_cols().is_compatible(*Gc.space_rows()->sub_space(equ_undecomp())) == true</tt>
   * <li> [<tt>Uy != NULL</tt>] <tt>Uy.space_rows().is_compatible(*space_range()) == true</tt>
   * <li> [<tt>Uy != NULL</tt>] <tt>Uy == Gc(:,equ_undecomp())'*Y</tt>
   * <li> [<tt>mat_rel == MATRICES_INDEP_IMPS</tt>] The behaviors of all of the participating output matrix
   *      objects must not be altered by changes to other matrix objects.
   * <li> [<tt>mat_rel == MATRICES_ALLOW_DEP_IMPS</tt>] The behaviors of all of the participating output matrix
   *      objects may change when other matrix objects or \c this is altered.
   * </ul>
   *
   * Note that this method requires that all of the output matrix objects \c Z, \c Y, \c R,
   * \c Uz, \c Uy, \c Vz and \c Vy must be independent of \c this and of each other.
   * For example, the behavior of \c R must be be altered if \c this is destroyed or if
   * \c Z is modified of destroyed.  This requirment constrains the implementations somewhat
   * but makes things much easier for the client and gives the client much more power.
   */
  virtual void update_decomp(
    std::ostream          *out
    ,EOutputLevel         olevel
    ,ERunTests            test_what
    ,const MatrixOp       &Gc
    ,MatrixOp             *Z
    ,MatrixOp             *Y
    ,MatrixOpNonsing      *R
    ,MatrixOp             *Uz
    ,MatrixOp             *Vy
    ,EMatRelations        mat_rel = MATRICES_INDEP_IMPS
    ) const = 0;
  
  /** \brief Print the sub-algorithm by which the decomposition is formed
   *
   */
  virtual void print_update_decomp(
    std::ostream& out, const std::string& leading_str ) const = 0;

  //@}
  
};	// end class DecompositionSystem

}	// end namespace ConstrainedOptPack

#endif // DECOMPOSITION_SYSTEM_H