MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_MatrixComposite.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 MATRIX_COMPOSITE_STD_H
43 #define MATRIX_COMPOSITE_STD_H
44 
45 #include <deque>
46 
50 #include "Teuchos_RCP.hpp"
51 #include "ReleaseResource.hpp"
52 
53 namespace AbstractLinAlgPack {
54 
126 class MatrixComposite : public MatrixOp {
127 public:
128 
129  // ///////////////////////////////////
130  // Public types
131 
133  typedef Teuchos::RCP<
135 
140  struct SubVectorEntry {
145  size_type r_l, size_type c_l, value_type beta
146  ,const Range1D& rng_G
147  ,const GPMS_ptr_t& G, const release_resource_ptr_t& G_release, BLAS_Cpp::Transp G_trans
148  ,const Vector* v
149  ,const release_resource_ptr_t& v_release, BLAS_Cpp::Transp v_trans
150  )
151  :r_l_(r_l),c_l_(c_l),beta_(beta),rng_G_(rng_G),G_(G),G_release_(G_release),G_trans_(G_trans)
152  ,v_(v),v_release_(v_release),v_trans_(v_trans)
153  {}
156  {
157  return
158  r_l_==v.r_l_ && c_l_==v.c_l_ && beta_==v.beta_
159  && rng_G_==v.rng_G_ && G_.get()==v.G_.get() && G_release_.get()==v.G_release_.get() && G_trans_==v.G_trans_
160  && v_==v.v_ && v_release_.get()==v.v_release_.get() && v_trans_==v.v_trans_;
161  }
164  c_l_;
165 
167 
169 
171 
175 
176  const Vector *v_;
177 
181  }; // end struct SubVectorEntry
182 
184  typedef std::deque<SubVectorEntry> vector_list_t;
185 
190  struct SubMatrixEntry {
195  size_type r_l, size_type r_u, size_type c_l, size_type c_u, value_type alpha
196  ,const Range1D& rng_P
197  ,const GPMS_ptr_t& P, const release_resource_ptr_t& P_release, BLAS_Cpp::Transp P_trans
198  ,const MatrixOp* A, const release_resource_ptr_t& A_release, BLAS_Cpp::Transp A_trans
199  ,const Range1D& rng_Q
200  ,const GPMS_ptr_t& Q, const release_resource_ptr_t& Q_release, BLAS_Cpp::Transp Q_trans
201  )
202  :r_l_(r_l),r_u_(r_u),c_l_(c_l),c_u_(c_u),alpha_(alpha),rng_P_(rng_P),P_(P),P_release_(P_release),P_trans_(P_trans)
203  ,A_(A),A_release_(A_release),A_trans_(A_trans),rng_Q_(rng_Q),Q_(Q),Q_release_(Q_release),Q_trans_(Q_trans)
204  {}
207  {
208  return
209  r_l_==m.r_l_ && r_u_==m.r_u_ && c_l_==m.c_l_ && c_u_==m.c_u_ && alpha_==m.alpha_
210  && rng_P_==m.rng_P_ && P_.get()==m.P_.get() && P_release_.get()==m.P_release_.get() && P_trans_==m.P_trans_
211  && A_==m.A_ && A_release_.get()==m.A_release_.get() && A_trans_==m.A_trans_
212  && rng_Q_==m.rng_Q_ && Q_.get()==m.Q_.get() && Q_release_.get()==m.Q_release_.get() && Q_trans_==m.Q_trans_;
213  }
219  Range1D rng_P_; // rng_P_.size() > 0 => P_ is ignored, rng_P_.full_range() => all rows op(A)
227  const MatrixOp *A_;
233  Range1D rng_Q_; // rng_Q_.size() > 0 => Q_ is ignored, rng_Q_.full_range() => all columns op(A)
240  }; // end struct SubMatrixEntry
241 
243  typedef std::deque<SubMatrixEntry> matrix_list_t;
244 
247 
253 
265  void reinitialize( size_type rows = 0, size_type cols = 0 );
266 
271  void add_vector(
272  size_type row_offset
273  ,size_type col_offset
274  ,value_type beta
275  ,const GenPermMatrixSlice *G
276  ,const release_resource_ptr_t &G_release
277  ,BLAS_Cpp::Transp G_trans
278  ,const Vector *v
279  ,const release_resource_ptr_t &v_release
280  ,BLAS_Cpp::Transp v_trans
281  );
282 
287  void add_vector(
288  size_type row_offset
289  ,size_type col_offset
290  ,value_type beta
291  ,const Range1D &rng_G
292  ,const Vector *v
293  ,const release_resource_ptr_t &v_release
294  ,BLAS_Cpp::Transp v_trans
295  );
296 
301  void add_vector(
302  size_type row_offset
303  ,size_type col_offset
304  ,value_type beta
305  ,const Vector *v
306  ,const release_resource_ptr_t &v_release
307  ,BLAS_Cpp::Transp v_trans
308  );
309 
317  void remove_vector( vector_list_t::iterator itr );
318 
323  void add_matrix(
324  size_type row_offset
325  ,size_type col_offset
326  ,value_type alpha
327  ,const GenPermMatrixSlice *P
328  ,const release_resource_ptr_t &P_release
329  ,BLAS_Cpp::Transp P_trans
330  ,const MatrixOp *A
331  ,const release_resource_ptr_t &A_release
332  ,BLAS_Cpp::Transp A_trans
333  ,const GenPermMatrixSlice *Q
334  ,const release_resource_ptr_t &Q_release
335  ,BLAS_Cpp::Transp Q_trans
336  );
337 
342  void add_matrix(
343  size_type row_offset
344  ,size_type col_offset
345  ,value_type alpha
346  ,const Range1D &rng_P
347  ,const MatrixOp *A
348  ,const release_resource_ptr_t &A_release
349  ,BLAS_Cpp::Transp A_trans
350  ,const Range1D &rng_Q
351  );
352 
357  void add_matrix(
358  size_type row_offset
359  ,size_type col_offset
360  ,value_type alpha
361  ,const Range1D &rng_P
362  ,const MatrixOp *A
363  ,const release_resource_ptr_t &A_release
364  ,BLAS_Cpp::Transp A_trans
365  ,const GenPermMatrixSlice *Q
366  ,const release_resource_ptr_t &Q_release
367  ,BLAS_Cpp::Transp Q_trans
368  );
369 
374  void add_matrix(
375  size_type row_offset
376  ,size_type col_offset
377  ,value_type alpha
378  ,const GenPermMatrixSlice *P
379  ,const release_resource_ptr_t &P_release
380  ,BLAS_Cpp::Transp P_trans
381  ,const MatrixOp *A
382  ,const release_resource_ptr_t &A_release
383  ,BLAS_Cpp::Transp A_trans
384  ,const Range1D &rng_Q
385  );
386 
391  void add_matrix(
392  size_type row_offset
393  ,size_type col_offset
394  ,value_type alpha
395  ,const MatrixOp *A
396  ,const release_resource_ptr_t &A_release
397  ,BLAS_Cpp::Transp A_trans
398  );
399 
404  void add_matrix(
405  size_type row_offset
406  ,size_type col_offset
407  ,value_type alpha
408  ,const GenPermMatrixSlice *P
409  ,const release_resource_ptr_t &P_release
410  ,BLAS_Cpp::Transp P_trans
411  );
412 
420  void remove_matrix( matrix_list_t::iterator itr );
421 
452  virtual void finish_construction(
455  );
456 
458 
461 
463  int num_vectors() const;
465  vector_list_t::iterator vectors_begin();
467  vector_list_t::iterator vectors_end();
469  vector_list_t::const_iterator vectors_begin() const;
471  vector_list_t::const_iterator vectors_end() const;
473  int num_matrices() const;
475  matrix_list_t::iterator matrices_begin();
477  matrix_list_t::iterator matrices_end();
479  matrix_list_t::const_iterator matrices_begin() const;
481  matrix_list_t::const_iterator matrices_end() const;
482 
484 
487 
489  size_type rows() const;
491  size_type cols() const;
493  size_type nz() const;
494 
496 
499 
501  const VectorSpace& space_rows() const;
503  const VectorSpace& space_cols() const;
505  mat_ptr_t sub_view(const Range1D& row_rng, const Range1D& col_rng) const;
507  void Vp_StMtV(VectorMutable* vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
508  , const Vector& v_rhs2, value_type beta) const;
510  void Vp_StMtV(VectorMutable* vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
511  , const SpVectorSlice& sv_rhs2, value_type beta) const;
513  void Vp_StPtMtV(VectorMutable* vs_lhs, value_type alpha
514  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
515  , BLAS_Cpp::Transp M_rhs2_trans
516  , const Vector& v_rhs3, value_type beta) const;
518  void Vp_StPtMtV(VectorMutable* vs_lhs, value_type alpha
519  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
520  , BLAS_Cpp::Transp M_rhs2_trans
521  , const SpVectorSlice& sv_rhs3, value_type beta) const;
522 
524 
525 private:
526 
527  // ///////////////////////////////
528  // private data members
529 
532 #ifdef DOXYGEN_COMPILE
537 #else
538  matrix_list_t matrix_list_;
539  vector_list_t vector_list_;
540  VectorSpace::space_ptr_t space_cols_;
541  VectorSpace::space_ptr_t space_rows_;
542 #endif
543 
544  // ///////////////////////////////
545  // private member functions
546 
547  void assert_fully_constructed() const;
548 
549 }; // end class MatrixComposite
550 
551 } // end namespace AbstractLinAlgPack
552 
553 #endif // MATRIX_COMPOSITE_STD_H
void remove_vector(vector_list_t::iterator itr)
Remove a sub-vector.
std::deque< SubMatrixEntry > matrix_list_t
Warning! This could be changed to some other STL container!
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
void add_matrix(size_type row_offset, size_type col_offset, value_type alpha, const GenPermMatrixSlice *P, const release_resource_ptr_t &P_release, BLAS_Cpp::Transp P_trans, const MatrixOp *A, const release_resource_ptr_t &A_release, BLAS_Cpp::Transp A_trans, const GenPermMatrixSlice *Q, const release_resource_ptr_t &Q_release, BLAS_Cpp::Transp Q_trans)
Add a sub-matrix alpha*op(P)*op(A)*op(Q).
mat_ptr_t sub_view(const Range1D &row_rng, const Range1D &col_rng) const
T * get() const
Teuchos::RCP< MemMngPack::ReleaseResource > release_resource_ptr_t
size_type c_l_
column of first element of vector in composite matrix.
virtual void finish_construction(const VectorSpace::space_ptr_t &space_cols, const VectorSpace::space_ptr_t &space_rows)
Call to finish the construction process.
std::deque< SubVectorEntry > vector_list_t
Warning! This could be changed to some other STL container!
. One-based subregion index range class.
Abstract interface for objects that represent a space for mutable coordinate vectors.
Abstract interface for releasing an object when it is not needed anymore {abstract}.
void Vp_StPtMtV(VectorMutable *vs_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const Vector &v_rhs3, value_type beta) const
const VectorSpace & space_cols() const
Base class for all matrices that support basic matrix operations.
MatrixComposite(size_type rows=0, size_type cols=0)
Construct.
void reinitialize(size_type rows=0, size_type cols=0)
Initialize a sized (on unsized) zero matrix to start with.
Abstract interface for mutable coordinate vectors {abstract}.
void remove_matrix(matrix_list_t::iterator itr)
Remove a sub-matrix.
SubMatrixEntry(size_type r_l, size_type r_u, size_type c_l, size_type c_u, value_type alpha, const Range1D &rng_P, const GPMS_ptr_t &P, const release_resource_ptr_t &P_release, BLAS_Cpp::Transp P_trans, const MatrixOp *A, const release_resource_ptr_t &A_release, BLAS_Cpp::Transp A_trans, const Range1D &rng_Q, const GPMS_ptr_t &Q, const release_resource_ptr_t &Q_release, BLAS_Cpp::Transp Q_trans)
Matrix class for matrices composed out of a set of other matrices and vectors.
SubVectorEntry(size_type r_l, size_type c_l, value_type beta, const Range1D &rng_G, const GPMS_ptr_t &G, const release_resource_ptr_t &G_release, BLAS_Cpp::Transp G_trans, const Vector *v, const release_resource_ptr_t &v_release, BLAS_Cpp::Transp v_trans)
Transp
TRANS.
void add_vector(size_type row_offset, size_type col_offset, value_type beta, const GenPermMatrixSlice *G, const release_resource_ptr_t &G_release, BLAS_Cpp::Transp G_trans, const Vector *v, const release_resource_ptr_t &v_release, BLAS_Cpp::Transp v_trans)
Add a sub-vector beta*op(op(G)*v).
const VectorSpace & space_rows() const
Concrete matrix type to represent general permutation (mapping) matrices.
void Vp_StMtV(VectorMutable *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const